Class axis_z_mod
In: shared/axis_z.f90

begin

Dependency

Methods

Included Modules

type_mod axis_type_mod nmlfile_mod grid_3d_mod constants_mod dc_types dc_url dc_string dc_trace dc_message gt4_history

Public Instance methods

Subroutine :

Dependency

[Source]

  subroutine axis_z_end
  !==== Dependency
    use dc_trace,  only: DbgMessage, BeginSub, EndSub
    use gt4_history, only: HistoryAxisClear
                                                                 !=end
    implicit none
    character(STRING), parameter:: subname = "axis_z_end"
  continue
    !----------------------------------------------------------------
    !   Check Initialization
    !----------------------------------------------------------------
    call BeginSub(subname)
    if ( .not. axis_z_initialized) then
       call EndSub( subname, 'axis_z_init was not called', c1=trim(subname) )
       return
    else
       axis_z_initialized = .false.
       axis_z_data_from_sigmahalf = .false.
       axis_z_data_from_manual   = .false.
       axis_z_data_from_netcdf   = ''
       axis_z_half_data_from_manual   = .false.
       axis_z_half_data_from_netcdf   = ''

       !--------------------------------------------------------------
       !   Initialize z_Dim
       !--------------------------------------------------------------
       z_Dim%stored   = .false.
       call HistoryAxisClear(z_Dim % axisinfo)
!!$       z_Dim%axisinfo%name     = ''
!!$       z_Dim%axisinfo%length   = 0
!!$       z_Dim%axisinfo%longname = ''
!!$       z_Dim%axisinfo%units    = ''
!!$       z_Dim%axisinfo%xtype    = ''
!!$       if ( associated(z_Dim%a_Dim) ) then
!!$          deallocate( z_Dim%a_Dim )
!!$       endif
!!$
!!$       if ( associated(z_Dim%attrs) ) then
!!$          deallocate( z_Dim%attrs )
!!$       end if

       !--------------------------------------------------------------
       !   Initialize z_DimHalf
       !--------------------------------------------------------------
       z_DimHalf%stored   = .false.
       call HistoryAxisClear(z_DimHalf % axisinfo)
!!$       z_DimHalf%axisinfo%name     = ''
!!$       z_DimHalf%axisinfo%length   = 0
!!$       z_DimHalf%axisinfo%longname = ''
!!$       z_DimHalf%axisinfo%units    = ''
!!$       z_DimHalf%axisinfo%xtype    = ''
!!$       if ( associated(z_DimHalf%a_Dim) ) then
!!$          deallocate( z_DimHalf%a_Dim )
!!$       endif
!!$
!!$       if ( associated(z_DimHalf%attrs) ) then
!!$          deallocate( z_DimHalf%attrs )
!!$       end if

    endif

    call EndSub( subname )
  end subroutine axis_z_end
Subroutine :
DimHalf :type(AXISINFO), intent(inout)
: 次元情報を包括する変数 =end

[Source]

  subroutine axis_z_half_manual(DimHalf)
    use axis_type_mod, only: axis_type_copy
    use dc_trace,   only: DbgMessage, BeginSub, EndSub
                                                                 !=end
    implicit none
                                                                 !=begin
    !==== In/Out
    !
    type(AXISINFO), intent(inout) :: DimHalf        ! 次元情報を包括する変数
                                                                 !=end
    character(STRING), parameter:: subname = "axis_z_half_manual"
  continue
    !-----------------------------------------------------------------
    !   Check Initialization
    !-----------------------------------------------------------------
    call BeginSub( subname )
    if (.not. axis_z_initialized) then
       call EndSub( subname, 'Call axis_z_init before call %c', c1=trim(subname) )
       return
    endif

    !----------------------------------------------------------------
    !   decision が 'manual' でない場合は停止して返す。
    !----------------------------------------------------------------
    if (.not. axis_z_half_data_from_manual) then
       call EndSub( subname, 'In axis_z_half_nml, configurated to not generate from manual')
       return
    endif

    !----------------------------------------------------------------
    !   intent(out) の引数にデータを渡す。
    !----------------------------------------------------------------
    call axis_type_copy(z_DimHalf, DimHalf)
    call EndSub(subname)
  end subroutine axis_z_half_manual
Subroutine :
DimHalf :type(AXISINFO), intent(inout)
: 次元情報を包括する変数 =end

Dependency

[Source]

  subroutine axis_z_half_netcdf(DimHalf)
  !==== Dependency
    use type_mod,      only: STRING
    use axis_type_mod, only: axis_type_copy
    use gt4_history,only: HistoryGet
    use dc_url  ,   only: UrlSplit
    use dc_trace,   only: DbgMessage, BeginSub, EndSub
                                                                 !=end
    implicit none
                                                                 !=begin
    !==== In/Out
    !
    type(AXISINFO), intent(inout) :: DimHalf        ! 次元情報を包括する変数
                                                                 !=end

    character(STRING)           :: file
    character(STRING)           :: varname

    character(STRING), parameter:: subname = "axis_z_half_netcdf"
  continue
    !-----------------------------------------------------------------
    !   Check Initialization
    !-----------------------------------------------------------------
    call BeginSub( subname )
    if (.not. axis_z_initialized) then
       call EndSub( subname, 'Call axis_z_init before call %c', c1=trim(subname) )
       return
    endif

    !----------------------------------------------------------------
    !   decision が '***@lon' などでない場合は停止して返す。
    !----------------------------------------------------------------
    if (axis_z_half_data_from_netcdf == '') then
       call EndSub( subname, 'In axis_z_half_nml, configurated to not generate from netcdf data')
       return
    endif

    !----------------------------------------------------------------
    !   HistoryGet によるデータ入力
    !----------------------------------------------------------------
    call UrlSplit(axis_z_half_data_from_netcdf, file=file, var=varname)
    call HistoryGet(file, varname, z_DimHalf%a_Dim)
    z_DimHalf%stored   = .true.

    !----------------------------------------------------------------
    !   intent(out) の引数にデータを渡す。
    !----------------------------------------------------------------
    call axis_type_copy(z_DimHalf, DimHalf)
    call EndSub(subname)
  end subroutine axis_z_half_netcdf
Subroutine :

Dependency

This procedure input/output NAMELIST#axis_z_nml, NAMELIST#axis_z_half_nml, NAMELIST#axis_z_attr_nml, NAMELIST#axis_z_half_attr_nml .

[Source]

  subroutine axis_z_init
  !==== Dependency
    use type_mod   ,only: STRING, TOKEN, INTKIND, REKIND, DBKIND, NMLARRAY
    use nmlfile_mod,only: nmlfile_init, nmlfile_open, nmlfile_close
    use grid_3d_mod,only: km, grid_3d_init
!    use axis_type_mod, only : axis_attrs_copy, axis_attrs_init
    use constants_mod,only: constants_init
!    use gt4_history, only: GT_HISTORY_ATTR
    use dc_types,  only: GT_TOKEN => TOKEN, GT_STRING => STRING
    use dc_url     ,only: GT_ATMARK, GT_QUESTION
    use dc_string  ,only: toChar
    use dc_trace   ,only: DbgMessage, BeginSub, EndSub
    use dc_message ,only: MessageNotify
    use gt4_history, only: HistoryAxisCreate, HistoryAxisAddAttr
    use dc_string, only: StrHead
                                                                 !=end
    implicit none

    !-------------------------------------------------------------------
    !   変数定義
    !-------------------------------------------------------------------
                                                                 !=begin
    integer(INTKIND)  :: z_Dim_length, z_DimHalf_length
    !
    !==== NAMELIST axis_z_nml
    !
    !axis_z_nml には、Z 軸の次元変数に関する情報を与える。
    !値を与えないものに関しては以下のデフォルトの値が用いられる。
    !
    !変数 decision には Z 軸のデータをどのように与えるかを指定する。
    !
    !* (({ 'manual' }))
    !  * Data 配列に格納したデータをそのまま Z 軸として与える。
    !
    !* (({ 'sigmahalf' }))
    !  * Z 軸が整数σレベルであると想定し、 axis_z_half_nml で
    !    与えられる半整数σレベルの値から、自動的に
    !    Z 軸のデータを決める。
    !    よって、axis_z_half_nml に有効なデータが与えられない
    !    場合にはデータは定まらない。
    !
    !* gtool4 変数 (例えば (({ 'foo.nc@lon' })) など)
    !  * 該当する変数から Z 軸のデータを取得する。
    !
    !* 上記以外
    !  * 現在は (({ 'sigmahalf' })) の場合と同様に設定される
    !
    !変数 length には、 ((<grid_3d_mod>)) の公開要素 ((< km >)) と同じ
    !値を与えなければならない。
    !
    character(TOKEN)  :: name     = 'sigma'       ! 次元変数名
    integer(INTKIND)  :: length   = 12            ! 次元長 (配列サイズ)
    character(STRING) :: longname = 'sigma at full level' ! 次元変数の記述的名称
    character(STRING) :: units    = 'sigma_level' ! 次元変数の単位
    character(TOKEN)  :: xtype    = 'float'       ! 次元変数の型
    character(STRING) :: decision = 'sigmahalf'   ! 次元データの取得方法
    real(REKIND)      :: Data(NMLARRAY)   = 0.0! 次元データ入力用

    namelist /axis_z_nml/ name         , length       , longname     , units        , xtype        , decision     , Data              ! 次元データ


    !==== NAMELIST axis_z_half_nml
    !
    !axiz_z_half_nml には、Z 軸の半整数レベルの次元変数に関する
    !情報を与える。
    !値を与えないものに関しては以下のデフォルトの値が用いられる。
    !
    !変数 decision には Z 軸のデータをどのように与えるかを指定する。
    !
    !* (({ 'manual' }))
    !  * Data 配列に格納したデータをそのまま 
    !    Z 軸半整数レベルとして与える。
    !
    !* gtool4 変数 (例えば (({ 'foo.nc@lon' })) など)
    !  * 該当する変数から Z 軸のデータを取得する。
    !
    !* 上記以外
    !  * 現在は (({ 'manual' })) の場合と同様に設定される
    !
    !変数 length には、 ((<grid_3d_mod>)) の公開要素 ((< km >)) に
    !プラス 1 した値を与えなければならない。
    !
!    name        = 'sigmahalf'   ! 次元変数名
!    length      = 13            ! 次元長 (配列サイズ)
!    longname    = 'sigma at half level' ! 次元変数の記述的名称
!    units       = 'sigma_level' ! 次元変数の単位
!    xtype       = 'float'       ! 次元変数の型
!    decision    = 'manual'      ! 次元データの取得方法
!    Data(1:13)  = (/1, 0.99, 0.97, 0.93, 0.85, 0.75, 0.63, 0.5, &
!         &             0.36, 0.22, 0.1, 0.05, 0/)  ! 次元データ入力用

    namelist /axis_z_half_nml/ name         , length       , longname     , units        , xtype        , decision     , Data              ! 次元データ
                                                                 !=end
                                                                 !=begin
    !==== NAMELIST axis_z_attr_nml
    !
    !Z 軸の次元変数の属性に関する情報を与える。
    !NAMELIST に複数の axis_z_attr_nml を用意しておく事で
    !複数の情報を与える事が可能である。
    !与えない場合には属性情報は付加されない。
    !
    !attrtype には与える属性値の種類を設定する。
    !((<URL:http://www.gfd-dennou.org/arch/gtool4/gt4f90io-current/doc/gt_history.htm#derived_gthistoryattr>))
    !を参照せよ。なお、arraysize に 1 以上の値を設定すると、
    !配列データが優先されて属性値に設定される。
    !
    character(GT_TOKEN)  :: attrname = '' ! 属性名
    character(GT_TOKEN)  :: attrtype = '' ! 属性値の型
    character(GT_STRING) :: cvalue   = '' ! 属性の値 (文字)
    integer(INTKIND)     :: ivalue   = 0      ! 属性の値 (整数)
    real(REKIND)         :: rvalue   = 0.0    ! 属性の値 (単精度実数)
    real(DBKIND)         :: dvalue   = 0.0d0  ! 属性の値 (倍精度実数)
    logical              :: lvalue   = .false.! 属性の値 (論理)
    integer(INTKIND)     :: arraysize= 0      ! 配列のサイズ
    integer(INTKIND) :: iarray(NMLARRAY)  = 0    ! 属性の値 (整数)
    real(REKIND)     :: rarray(NMLARRAY)  = 0.0  ! 属性の値 (単精度実数)
    real(DBKIND)     :: darray(NMLARRAY)  = 0.0d0! 属性の値 (倍精度実数)

    namelist /axis_z_attr_nml/ attrname     , attrtype     , cvalue       , ivalue       , rvalue       , dvalue       , lvalue       , arraysize    , iarray       , rarray       , darray            ! 属性の値 (倍精度実数)

    !==== NAMELIST axis_z_half_attr_nml
    !
    !Z 軸の半整数レベルの次元変数の属性に関する情報を与える。
    !NAMELIST に複数の axis_z_half_attr_nml を用意しておく事で
    !複数の情報を与える事が可能である。
    !与えない場合には属性情報は付加されない。
    !
    !attrtype には与える属性値の種類を設定する。
    !((<URL:http://www.gfd-dennou.org/arch/gtool4/gt4f90io-current/doc/gt_history.htm#derived_gthistoryattr>))
    !を参照せよ。なお、arraysize に 1 以上の値を設定すると、
    !配列データが優先されて属性値に設定される。
    !
    namelist /axis_z_half_attr_nml/ attrname     , attrtype     , cvalue       , ivalue       , rvalue       , dvalue       , lvalue       , arraysize    , iarray       , rarray       , darray            ! 属性の値 (倍精度実数)
                                                                 !=end


    !-----------------------------------------------------------------
    !   変数情報の一時格納変数
    !-----------------------------------------------------------------
!    type(GT_HISTORY_ATTR), allocatable :: attrs_tmp(:)

    !-----------------------------------------------------------------
    !   汎用変数
    !-----------------------------------------------------------------
    integer(INTKIND)            :: i, k
    integer(INTKIND)            :: nmlstat, nmlunit
    logical                     :: nmlreadable, next
    character(TOKEN)            :: position

    character(STRING), parameter:: subname = "axis_z_init"
  continue
    !----------------------------------------------------------------
    !   Check Initialization
    !----------------------------------------------------------------
    call BeginSub(subname)
    if (axis_z_initialized) then
       call EndSub( subname, '%c is already called', c1=trim(subname) )
       return
    else
       axis_z_initialized = .true.
    endif

    !----------------------------------------------------------------
    !   Version identifier
    !----------------------------------------------------------------
    call DbgMessage('%c :: %c', c1=trim(version), c2=trim(tagname))

    !----------------------------------------------------------------
    !   read axis_z_nml
    !----------------------------------------------------------------
    call nmlfile_init
    call nmlfile_open(nmlunit, nmlreadable)
    if (nmlreadable) then
       read(nmlunit, nml=axis_z_nml, iostat=nmlstat)
       call DbgMessage('Stat of NAMELIST axis_z_nml Input is <%d>', i=(/nmlstat/))
       write(0, nml=axis_z_nml)
    else
       call DbgMessage('Not Read NAMELIST axis_z_nml')
       call MessageNotify('W', subname, 'Can not Read NAMELIST axis_z_nml. Force Use Default Value.')
    end if
    call nmlfile_close

    z_Dim%stored   = .false.
    !  次元変数の情報を構造型変数 z_Dim への代入
    z_Dim_length = length
    call HistoryAxisCreate(z_Dim % axisinfo, name, z_Dim_length, longname, units, xtype)
!!$    z_Dim%axisinfo%name     = name
!!$    z_Dim%axisinfo%length   = length
!!$    z_Dim%axisinfo%longname = longname
!!$    z_Dim%axisinfo%units    = units
!!$    z_Dim%axisinfo%xtype    = xtype
    allocate( z_Dim%a_Dim(z_Dim_length) )
    call DbgMessage('z_Dim-length=<%d>', i=(/z_Dim_length/))


    !  次元変数の情報を構造型変数 z_Dim への代入
    select case(decision)

    ! manual: NAMELIST の Data で入力
    case('manual')
       z_Dim%a_Dim(:)                       = 0
       z_Dim%a_Dim(1:z_Dim_length) = Data(1:z_Dim_length)
       z_Dim%stored                         = .true.
       axis_z_data_from_manual              = .true.

    ! sigmahalf: 半整数レベルより作成
    case('sigmahalf')
       axis_z_data_from_sigmahalf = .true.
       z_Dim%stored               = .false.

    ! foo.nc@lon: foo.nc ファイルの lon 変数から取得
    ! その他    : sigmahalf と同じに
    case default
       ! 文字の中に '@' か '?' が含まれる場合は gtool4 変数として
       ! 認識し、その nc ファイルから変数情報をコピーする。
       if (   index(decision, GT_ATMARK)   > 0 .or. index(decision, GT_QUESTION) > 0) then
          axis_z_data_from_netcdf  = decision
          z_Dim%stored             = .false.
       ! それ以外は 'sigmahalf' と同じように処理
       else
          axis_z_data_from_sigmahalf = .true.
          z_Dim%stored         = .false.
       endif
    end select

    !----------------------------------------------------------------
    !   read axis_z_attr_nml
    !----------------------------------------------------------------
    call nmlfile_init
    call nmlfile_open(nmlunit, nmlreadable)

    if (.not. nmlreadable) then
       call DbgMessage('Not Read NAMELIST axis_z_attr_nml')
       call MessageNotify('W', subname, 'Can not Read NAMELIST axis_z_attr_nml.')
    else

       i = 0
       next = .false.
       axis_z_attr_nml_input : do
          i = i + 1
          call DbgMessage('NAMELIST axis_z_attr_nml Input, <%d> time', i=(/i/))
          ! 初期化
          attrname  = ''     ! 属性名
          attrtype  = ''     ! 属性値の型
          cvalue    = ''     ! 属性の値 (文字)
          ivalue    = 0      ! 属性の値 (整数)
          rvalue    = 0.0    ! 属性の値 (単精度実数)
          dvalue    = 0.0d0  ! 属性の値 (倍精度実数)
          lvalue    = .false.! 属性の値 (論理)
          arraysize = 0      ! 配列のサイズ
          iarray(:) = 0      ! 属性の値 (整数)
          rarray(:) = 0.0    ! 属性の値 (単精度実数)
          darray(:) = 0.0d0  ! 属性の値 (倍精度実数)

          ! read nml
          read(nmlunit, nml=axis_z_attr_nml, iostat=nmlstat)
          call DbgMessage('Stat of NAMELIST axis_z_attr_nml Input is <%d>', i=(/nmlstat/))
          write(0, nml=axis_z_attr_nml)

          ! Inquire access position
          inquire(nmlunit, position=position)
          if ( trim(position) /= 'APPEND' .and. nmlstat == 0 ) then
             next = .true.
          else
             next = .false.
          endif

          ! 有効でない値を含むものに関しては無視。
          if (attrname == '') then
             call DbgMessage('attrname is blank. so this axis_z_attr_nml is ignored.')
             if (next) cycle
             if (.not. next) exit
          elseif (attrtype == '') then
             call DbgMessage('attrtype is blank. so this axis_z_attr_nml is ignored.')
             if (next) cycle
             if (.not. next) exit
          endif

          !-----------------------------------------------------------
          ! z_Dim%attrs への格納
          !-----------------------------------------------------------
          ! attrs(:) の拡張
!!$          if ( .not. associated(z_Dim%attrs) ) then
!!$             allocate( z_Dim%attrs(1) )
!!$             k = 1
!!$          else
!!$             k = size( z_Dim%attrs ) + 1
!!$             ! 配列データの領域確保
!!$             allocate( attrs_tmp(k-1) )
!!$             call axis_attrs_copy(from=z_Dim%attrs(1:k-1), to=attrs_tmp(1:k-1))
!!$             deallocate( z_Dim%attrs )
!!$             allocate( z_Dim%attrs(k) )
!!$             call axis_attrs_copy(from=attrs_tmp(1:k-1), to=z_Dim%attrs(1:k-1))
!!$             deallocate( attrs_tmp )
!!$          endif
!!$
!!$          if (arraysize > 0) then
!!$             call axis_attrs_init(z_Dim%attrs(k))
!!$
!!$             deallocate(  z_Dim%attrs(k)%iarray  )
!!$             deallocate(  z_Dim%attrs(k)%rarray  )
!!$             deallocate(  z_Dim%attrs(k)%darray  )
!!$
!!$             allocate(  z_Dim%attrs(k)%iarray( arraysize )  )
!!$             allocate(  z_Dim%attrs(k)%rarray( arraysize )  )
!!$             allocate(  z_Dim%attrs(k)%darray( arraysize )  )
!!$
!!$             z_Dim%attrs(k)%array = .true.
!!$
!!$          else
!!$             call axis_attrs_init(z_Dim%attrs(k))
!!$          endif

          if (StrHead(attrtype, 'char')) then
            call HistoryAxisAddAttr(z_Dim % axisinfo, attrname, cvalue)
          elseif (StrHead(attrtype, 'logical')) then
            call HistoryAxisAddAttr(z_Dim % axisinfo, attrname, lvalue)
          elseif (StrHead(attrtype, 'int')) then
            if (arraysize > 0) then
              call HistoryAxisAddAttr(z_Dim % axisinfo, attrname, iarray)
            else
              call HistoryAxisAddAttr(z_Dim % axisinfo, attrname, ivalue)
            end if
          elseif (StrHead(attrtype, 'real')) then
            if (arraysize > 0) then
              call HistoryAxisAddAttr(z_Dim % axisinfo, attrname, rarray)
            else
              call HistoryAxisAddAttr(z_Dim % axisinfo, attrname, rvalue)
            end if
          elseif (StrHead(attrtype, 'double')) then
            if (arraysize > 0) then
              call HistoryAxisAddAttr(z_Dim % axisinfo, attrname, darray)
            else
              call HistoryAxisAddAttr(z_Dim % axisinfo, attrname, dvalue)
            end if
          end if

!!$          z_Dim%attrs(k)%attrname  = attrname
!!$          z_Dim%attrs(k)%attrtype  = attrtype
!!$          z_Dim%attrs(k)%cvalue    = cvalue
!!$          z_Dim%attrs(k)%ivalue    = ivalue
!!$          z_Dim%attrs(k)%rvalue    = rvalue
!!$          z_Dim%attrs(k)%dvalue    = dvalue
!!$          z_Dim%attrs(k)%lvalue    = lvalue
!!$          z_Dim%attrs(k)%iarray(1:max(1,arraysize)) = iarray(1:max(1,arraysize))
!!$          z_Dim%attrs(k)%rarray(1:max(1,arraysize)) = rarray(1:max(1,arraysize))
!!$          z_Dim%attrs(k)%darray(1:max(1,arraysize)) = darray(1:max(1,arraysize))
!!$
!!$          call DbgMessage('z_Dim-attrs(%d) [attrname=<%c> '            // &
!!$               & 'attrtype=<%c> array=<%b> cvalue=<%c>  '         // &
!!$               & 'ivalue=<%d> rvalue=<%r> dvalue=<%f> '           // &
!!$               & 'iarray(1:%d)=<%d, ...> '                        // &
!!$               & 'rarray(1:%d)=<%r, ...> darray(1:%d)=<%f, ...>'   , &
!!$               & c1=trim( z_Dim%attrs(k)%attrname )                  , &
!!$               & c2=trim( z_Dim%attrs(k)%attrtype )                  , &
!!$               & c3=trim( z_Dim%attrs(k)%cvalue )                    , &
!!$               & i=(/ k, z_Dim%attrs(k)%ivalue    ,                    &
!!$               &      size(z_Dim%attrs(k)%iarray) ,                    &
!!$               &      z_Dim%attrs(k)%iarray       ,                    &
!!$               &      size(z_Dim%attrs(k)%rarray) ,                    &
!!$               &      size(z_Dim%attrs(k)%darray)                      &
!!$               &   /)                                                , &
!!$               & r=(/z_Dim%attrs(k)%rvalue, z_Dim%attrs(k)%rarray/)    , &
!!$               & d=(/z_Dim%attrs(k)%dvalue, z_Dim%attrs(k)%darray/)    , &
!!$               & l=(/z_Dim%attrs(k)%lvalue/)                      )

          if (.not. next) exit axis_z_attr_nml_input
          next      = .false.  ! 次回のための初期化

       enddo axis_z_attr_nml_input
    end if
    call nmlfile_close


    !----------------------------------------------------------------
    !   read axis_z_half_nml
    !----------------------------------------------------------------
    ! Initialize
    name        = 'sigmahalf'   ! 次元変数名
    length      = 13            ! 次元長 (配列サイズ)
    longname    = 'sigma at half level' ! 次元変数の記述的名称
    units       = 'sigma_level' ! 次元変数の単位
    xtype       = 'float'       ! 次元変数の型
    decision    = 'manual'      ! 次元データの取得方法
    Data(1:13)  = (/1.0, 0.99, 0.97, 0.93, 0.85, 0.75, 0.63, 0.5, 0.36, 0.22, 0.1, 0.05, 0.0/)  ! 次元データ入力用

    call nmlfile_init
    call nmlfile_open(nmlunit, nmlreadable)
    if (nmlreadable) then
       read(nmlunit, nml=axis_z_half_nml, iostat=nmlstat)
       call DbgMessage('Stat of NAMELIST axis_z_half_nml Input is <%d>', i=(/nmlstat/))
       write(0, nml=axis_z_half_nml)
    else
       call DbgMessage('Not Read NAMELIST axis_z_half_nml')
       call MessageNotify('W', subname, 'Can not Read NAMELIST axis_z_half_nml. Force Use Default Value.')
    end if
    call nmlfile_close

    z_DimHalf%stored   = .false.
    !  次元変数の情報を構造型変数 z_DimHalf への代入
    z_DimHalf_length = length
    call HistoryAxisCreate(z_DimHalf % axisinfo, name, z_DimHalf_length, longname, units, xtype)
!!$    z_DimHalf%axisinfo%name     = name
!!$    z_DimHalf%axisinfo%length   = length
!!$    z_DimHalf%axisinfo%longname = longname
!!$    z_DimHalf%axisinfo%units    = units
!!$    z_DimHalf%axisinfo%xtype    = xtype
    allocate( z_DimHalf%a_Dim(z_DimHalf_length) )
    call DbgMessage('z_DimHalf-length=<%d>', i=(/z_DimHalf_length/))


    !  次元変数の情報を構造型変数 z_DimHalf への代入
    select case(decision)

    ! manual: NAMELIST の Data で入力
    case('manual')
       z_DimHalf%a_Dim(:)           = 0
       z_DimHalf%a_Dim(1:z_DimHalf_length) = Data(1:z_DimHalf_length)
       z_DimHalf%stored             = .true.
       axis_z_half_data_from_manual = .true.

       call DbgMessage('z_DimHalf(:)=<%c>'                      , c1=trim( toChar(z_DimHalf%a_Dim(:)) )  )

    ! foo.nc@lon: foo.nc ファイルの lon 変数から取得
    ! その他    : sigmahalf と同じに
    case default
       ! 文字の中に '@' か '?' が含まれる場合は gtool4 変数として
       ! 認識し、その nc ファイルから変数情報をコピーする。
       if (   index(decision, GT_ATMARK)   > 0 .or. index(decision, GT_QUESTION) > 0) then
          axis_z_half_data_from_netcdf  = decision
          z_DimHalf%stored              = .false.
       ! それ以外は 'manual' と同じように処理
       else
          z_DimHalf%a_Dim(:)           = 0
          z_DimHalf%a_Dim(1:z_DimHalf_length) = Data(1:z_DimHalf_length)
          z_DimHalf%stored             = .true.
          axis_z_half_data_from_manual = .true.

          call DbgMessage('z_DimHalf(:)=<%c>'                      , c1=trim( toChar(z_DimHalf%a_Dim(:)) )  )
       endif
    end select

    !----------------------------------------------------------------
    !   read axis_z_half_attr_nml
    !----------------------------------------------------------------
    call nmlfile_init
    call nmlfile_open(nmlunit, nmlreadable)

    if (.not. nmlreadable) then
       call DbgMessage('Not Read NAMELIST axis_z_half_attr_nml')
       call MessageNotify('W', subname, 'Can not Read NAMELIST axis_z_half_attr_nml.')
    else

       i = 0
       next = .false.
       axis_z_half_attr_nml_input : do
          call DbgMessage('NAMELIST axis_z_half_attr_nml Input, <%d> time', i=(/i/))
          i = i + 1
          ! 初期化
          attrname  = ''     ! 属性名
          attrtype  = ''     ! 属性値の型
          cvalue    = ''     ! 属性の値 (文字)
          ivalue    = 0      ! 属性の値 (整数)
          rvalue    = 0.0    ! 属性の値 (単精度実数)
          dvalue    = 0.0d0  ! 属性の値 (倍精度実数)
          lvalue    = .false.! 属性の値 (論理)
          arraysize = 0      ! 配列のサイズ
          iarray(:) = 0      ! 属性の値 (整数)
          rarray(:) = 0.0    ! 属性の値 (単精度実数)
          darray(:) = 0.0d0  ! 属性の値 (倍精度実数)

          ! read nml
          read(nmlunit, nml=axis_z_half_attr_nml, iostat=nmlstat)
          call DbgMessage('Stat of NAMELIST axis_z_half_attr_nml Input is <%d>', i=(/nmlstat/))
          write(0, nml=axis_z_half_attr_nml)

          ! Inquire access position
          inquire(nmlunit, position=position)
          if ( trim(position) /= 'APPEND' .and. nmlstat == 0 ) then
             next = .true.
          else
             next = .false.
          endif

          ! 有効でない値を含むものに関しては無視。
          if (attrname == '') then
             call DbgMessage('attrname is blank. so this axis_z_half_attr_nml is ignored.')
             if (next) cycle
             if (.not. next) exit
          elseif (attrtype == '') then
             call DbgMessage('attrtype is blank. so this axis_z_half_attr_nml is ignored.')
             if (next) cycle
             if (.not. next) exit
          endif

          !-----------------------------------------------------------
          ! z_DimHalf%attrs への格納
          !-----------------------------------------------------------
          ! attrs(:) の拡張
!!$          if ( .not. associated(z_DimHalf%attrs) ) then
!!$             allocate( z_DimHalf%attrs(1) )
!!$             k = 1
!!$          else
!!$             k = size( z_DimHalf%attrs ) + 1
!!$             ! 配列データの領域確保
!!$             allocate( attrs_tmp(k-1) )
!!$             call axis_attrs_copy(from=z_DimHalf%attrs(1:k-1), to=attrs_tmp(1:k-1))
!!$             deallocate( z_DimHalf%attrs )
!!$             allocate( z_DimHalf%attrs(k) )
!!$             call axis_attrs_copy(from=attrs_tmp(1:k-1), to=z_DimHalf%attrs(1:k-1))
!!$             deallocate( attrs_tmp )
!!$          endif
!!$
!!$          if (arraysize > 0) then
!!$             call axis_attrs_init(z_DimHalf%attrs(k))
!!$
!!$             deallocate(  z_DimHalf%attrs(k)%iarray  )
!!$             deallocate(  z_DimHalf%attrs(k)%rarray  )
!!$             deallocate(  z_DimHalf%attrs(k)%darray  )
!!$
!!$             allocate(  z_DimHalf%attrs(k)%iarray( arraysize )  )
!!$             allocate(  z_DimHalf%attrs(k)%rarray( arraysize )  )
!!$             allocate(  z_DimHalf%attrs(k)%darray( arraysize )  )
!!$
!!$             z_DimHalf%attrs(k)%array = .true.
!!$
!!$          else
!!$             call axis_attrs_init(z_DimHalf%attrs(k))
!!$          endif

          if (StrHead(attrtype, 'char')) then
            call HistoryAxisAddAttr(z_DimHalf % axisinfo, attrname, cvalue)
          elseif (StrHead(attrtype, 'logical')) then
            call HistoryAxisAddAttr(z_DimHalf % axisinfo, attrname, lvalue)
          elseif (StrHead(attrtype, 'int')) then
            if (arraysize > 0) then
              call HistoryAxisAddAttr(z_DimHalf % axisinfo, attrname, iarray)
            else
              call HistoryAxisAddAttr(z_DimHalf % axisinfo, attrname, ivalue)
            end if
          elseif (StrHead(attrtype, 'real')) then
            if (arraysize > 0) then
              call HistoryAxisAddAttr(z_DimHalf % axisinfo, attrname, rarray)
            else
              call HistoryAxisAddAttr(z_DimHalf % axisinfo, attrname, rvalue)
            end if
          elseif (StrHead(attrtype, 'double')) then
            if (arraysize > 0) then
              call HistoryAxisAddAttr(z_DimHalf % axisinfo, attrname, darray)
            else
              call HistoryAxisAddAttr(z_DimHalf % axisinfo, attrname, dvalue)
            end if
          end if

!!$          z_DimHalf%attrs(k)%attrname  = attrname
!!$          z_DimHalf%attrs(k)%attrtype  = attrtype
!!$          z_DimHalf%attrs(k)%cvalue    = cvalue
!!$          z_DimHalf%attrs(k)%ivalue    = ivalue
!!$          z_DimHalf%attrs(k)%rvalue    = rvalue
!!$          z_DimHalf%attrs(k)%dvalue    = dvalue
!!$          z_DimHalf%attrs(k)%lvalue    = lvalue
!!$          z_DimHalf%attrs(k)%iarray(1:max(1,arraysize)) = iarray(1:max(1,arraysize))
!!$          z_DimHalf%attrs(k)%rarray(1:max(1,arraysize)) = rarray(1:max(1,arraysize))
!!$          z_DimHalf%attrs(k)%darray(1:max(1,arraysize)) = darray(1:max(1,arraysize))
!!$
!!$          call DbgMessage('z_DimHalf-attrs(%d) [attrname=<%c> '            // &
!!$               & 'attrtype=<%c> array=<%b> cvalue=<%c>  '         // &
!!$               & 'ivalue=<%d> rvalue=<%r> dvalue=<%f> '           // &
!!$               & 'iarray(1:%d)=<%d, ...> '                        // &
!!$               & 'rarray(1:%d)=<%r, ...> darray(1:%d)=<%f, ...>'   , &
!!$               & c1=trim( z_DimHalf%attrs(k)%attrname )                  , &
!!$               & c2=trim( z_DimHalf%attrs(k)%attrtype )                  , &
!!$               & c3=trim( z_DimHalf%attrs(k)%cvalue )                    , &
!!$               & i=(/ k, z_DimHalf%attrs(k)%ivalue    ,                    &
!!$               &      size(z_DimHalf%attrs(k)%iarray) ,                    &
!!$               &      z_DimHalf%attrs(k)%iarray       ,                    &
!!$               &      size(z_DimHalf%attrs(k)%rarray) ,                    &
!!$               &      size(z_DimHalf%attrs(k)%darray)                      &
!!$               &   /)                                                , &
!!$               & r=(/z_DimHalf%attrs(k)%rvalue, z_DimHalf%attrs(k)%rarray/)    , &
!!$               & d=(/z_DimHalf%attrs(k)%dvalue, z_DimHalf%attrs(k)%darray/)    , &
!!$               & l=(/z_DimHalf%attrs(k)%lvalue/)                      )

          if (.not. next) exit axis_z_half_attr_nml_input
          next      = .false.  ! 次回のための初期化

       enddo axis_z_half_attr_nml_input
    end if
    call nmlfile_close


    !----------------------------------------------------------------
    !   grid_3d_mod との整合性をチェック
    !----------------------------------------------------------------
    call grid_3d_init
    if (z_Dim_length /= km) then
       call MessageNotify('E', subname, message='axis length is inconsistent with km in grid_3d_mod')
    endif

    if (z_DimHalf_length /= z_Dim_length + 1) then
       call MessageNotify('E', subname, message='Sigma length should be SigmaHalf length - 1.')
    endif

    !-----------------------------------------------------------------
    !   constants モジュールの初期化
    !-----------------------------------------------------------------
    call constants_init

    !----------------------------------------------------------------
    !   例外処理
    !----------------------------------------------------------------
    if (length < 1) then
       call MessageNotify('E', subname, message='Invalid grid number.')
    endif

    call EndSub( subname )
  end subroutine axis_z_init
Subroutine :
Dim :type(AXISINFO), intent(inout)
: 次元情報を包括する変数 =end

[Source]

  subroutine axis_z_manual(Dim)
    use axis_type_mod, only: axis_type_copy
    use dc_trace,   only: DbgMessage, BeginSub, EndSub
                                                                 !=end
    implicit none
                                                                 !=begin
    !==== In/Out
    !
    type(AXISINFO), intent(inout) :: Dim        ! 次元情報を包括する変数
                                                                 !=end
    character(STRING), parameter:: subname = "axis_z_manual"
  continue
    !-----------------------------------------------------------------
    !   Check Initialization
    !-----------------------------------------------------------------
    call BeginSub( subname )
    if (.not. axis_z_initialized) then
       call EndSub( subname, 'Call axis_z_init before call %c', c1=trim(subname) )
       return
    endif

    !----------------------------------------------------------------
    !   decision が 'manual' でない場合は停止して返す。
    !----------------------------------------------------------------
    if (.not. axis_z_data_from_manual) then
       call EndSub( subname, 'In axis_z_nml, configurated to not generate from manual')
       return
    endif

    !----------------------------------------------------------------
    !   intent(out) の引数にデータを渡す。
    !----------------------------------------------------------------
    call axis_type_copy(z_Dim, Dim)
    call EndSub(subname)
  end subroutine axis_z_manual
Subroutine :
Dim :type(AXISINFO), intent(inout)
: 次元情報を包括する変数 =end

Dependency

[Source]

  subroutine axis_z_netcdf(Dim)
  !==== Dependency
    use type_mod,      only: STRING
    use axis_type_mod, only: axis_type_copy
    use gt4_history,only: HistoryGet
    use dc_url  ,   only: UrlSplit
    use dc_trace,   only: DbgMessage, BeginSub, EndSub
                                                                 !=end
    implicit none
                                                                 !=begin
    !==== In/Out
    !
    type(AXISINFO), intent(inout) :: Dim        ! 次元情報を包括する変数
                                                                 !=end

    character(STRING)           :: file
    character(STRING)           :: varname

    character(STRING), parameter:: subname = "axis_z_netcdf"
  continue
    !-----------------------------------------------------------------
    !   Check Initialization
    !-----------------------------------------------------------------
    call BeginSub( subname )
    if (.not. axis_z_initialized) then
       call EndSub( subname, 'Call axis_z_init before call %c', c1=trim(subname) )
       return
    endif

    !----------------------------------------------------------------
    !   decision が '***@lon' などでない場合は停止して返す。
    !----------------------------------------------------------------
    if (axis_z_data_from_netcdf == '') then
       call EndSub( subname, 'In axis_z_nml, configurated to not generate from netcdf data')
       return
    endif

    !----------------------------------------------------------------
    !   HistoryGet によるデータ入力
    !----------------------------------------------------------------
    call UrlSplit(axis_z_data_from_netcdf, file=file, var=varname)
    call HistoryGet(file, varname, z_Dim%a_Dim)
    z_Dim%stored   = .true.

    !----------------------------------------------------------------
    !   intent(out) の引数にデータを渡す。
    !----------------------------------------------------------------
    call axis_type_copy(z_Dim, Dim)
    call EndSub(subname)
  end subroutine axis_z_netcdf
Subroutine :
Dim :type(AXISINFO), intent(inout)
: 次元情報を包括する変数 =end
DimHalf :type(AXISINFO), intent(inout)
: 次元情報を包括する変数 =end

Dependency

[Source]

  subroutine axis_z_sigmahalf_manual(Dim, DimHalf)
  !==== Dependency
    use axis_type_mod, only: axis_type_copy
    use dc_trace,   only: DbgMessage, BeginSub, EndSub
                                                                 !=end
    implicit none
                                                                 !=begin
    !==== In/Out
    !
    type(AXISINFO), intent(inout) :: Dim, DimHalf  ! 次元情報を包括する変数
                                                                 !=end
    character(STRING), parameter:: subname = "axis_z_sigmahalf_manual"
  continue
    !-----------------------------------------------------------------
    !   Check Initialization
    !-----------------------------------------------------------------
    call BeginSub( subname )
    if (.not. axis_z_initialized) then
       call EndSub( subname, 'Call axis_z_init before call %c', c1=trim(subname) )
       return
    endif

    !----------------------------------------------------------------
    !   axis_z_nml の decision が 'sigmahalf' でなく、
    !   axis_z_sigmahalf_nml の decision が 'manual' でない場合は停止
    !----------------------------------------------------------------
    if (   .not. axis_z_data_from_sigmahalf .or. .not. axis_z_half_data_from_manual ) then
       call EndSub( subname, 'In axis_z_nml and in axis_z_sigmahalf_nml, '  // 'configurated to not generate from sigmahalf ' // '(which is generated from NAMELIST)')
       return
    endif

    !----------------------------------------------------------------
    !   DimHalf の値から Dim を作成。
    !----------------------------------------------------------------
    call gen_sigma_from_half(z_DimHalf, z_Dim)

    !----------------------------------------------------------------
    !   intent(out) の引数にデータを渡す。
    !----------------------------------------------------------------
    call axis_type_copy(z_Dim, Dim)
    call axis_type_copy(z_DimHalf, DimHalf)
    call EndSub(subname)
  end subroutine axis_z_sigmahalf_manual
Subroutine :
Dim :type(AXISINFO), intent(inout)
: 次元情報を包括する変数 =end
DimHalf :type(AXISINFO), intent(inout)
: 次元情報を包括する変数 =end

Dependency

[Source]

  subroutine axis_z_sigmahalf_netcdf(Dim, DimHalf)
  !==== Dependency
    use axis_type_mod, only: axis_type_copy
    use gt4_history,only: HistoryGet
    use dc_url  ,   only: UrlSplit
    use dc_trace,   only: DbgMessage, BeginSub, EndSub
                                                                 !=end
    implicit none
                                                                 !=begin
    !==== In/Out
    !
    type(AXISINFO), intent(inout) :: Dim, DimHalf  ! 次元情報を包括する変数
                                                                 !=end
    character(STRING)           :: file
    character(STRING)           :: varname

    character(STRING), parameter:: subname = "axis_z_sigmahalf_netcdf"
  continue
    !-----------------------------------------------------------------
    !   Check Initialization
    !-----------------------------------------------------------------
    call BeginSub( subname )
    if (.not. axis_z_initialized) then
       call EndSub( subname, 'Call axis_z_init before call %c', c1=trim(subname) )
       return
    endif

    !----------------------------------------------------------------
    !   axis_z_nml の decision が 'sigmahalf' でなく、
    !   axis_z_sigmahalf_nml の decision が '***@lon' でない時は停止
    !----------------------------------------------------------------
    if (   .not. axis_z_data_from_sigmahalf .or. .not. axis_z_half_data_from_netcdf /= '') then
       call EndSub( subname, 'In axis_z_nml and in axis_z_sigmahalf_nml, '  // 'configurated to not generate from sigmahalf ' // '(which is generated from netcdf data)')
       return
    endif

    !----------------------------------------------------------------
    !   HistoryGet によるデータ入力
    !----------------------------------------------------------------
    call UrlSplit(axis_z_half_data_from_netcdf, file=file, var=varname)
    call HistoryGet(file, varname, z_DimHalf%a_Dim)
    z_Dim%stored   = .true.

    !----------------------------------------------------------------
    !   DimHalf の値から Dim を作成。
    !----------------------------------------------------------------
    call gen_sigma_from_half(z_DimHalf, z_Dim)

    !----------------------------------------------------------------
    !   intent(out) の引数にデータを渡す。
    !----------------------------------------------------------------
    call axis_type_copy(z_Dim, Dim)
    call axis_type_copy(z_DimHalf, DimHalf)
    call EndSub(subname)
  end subroutine axis_z_sigmahalf_netcdf

[Validate]