gdncattrgetnum.f90

Path: gtdata/gtdata_netcdf/gdncattrgetnum.f90
Last Update: Tue Jun 22 23:13:47 +0900 2010

Get GD_NC_VARIABLES

Authors:Eizi TOYODA, Yasuhiro MORIKAWA
Version:$Id: gdncattrgetnum.f90,v 1.3 2010-06-22 14:13:47 morikawa Exp $
Tag Name:$Name: gtool5-20101228-1 $
Copyright:Copyright (C) GFD Dennou Club, 2000-2009. All rights reserved.
License:See COPYRIGHT

お客様向きではないけれど、情報落ちのないインターフェイスということで.…

stat < 0 :エラー、あるいはその属性は存在しなかった
stat = 0 … size(value):その属性を全部読み取った。サイズは stat 個
stat > size(value) :配列長不足のため属性が全部読み取れなかった。 サイズは stat 個必要
  • バグ:
    • 属性が文字型で STRING 文字を越える場合、GT_ECHARSHORT が返る

Required files

Methods

Included Modules

gtdata_netcdf_types gtdata_netcdf_internal netcdf_f77 dc_url gtdata_netcdf_generic dc_types dc_string

Public Instance methods

Subroutine :
var :type(GD_NC_VARIABLE), intent(in)
name :character(len = *), intent(in)
value(:) :real(DP), intent(out)
stat :integer, intent(out)
default :real(DP), intent(in), optional

[Source]

subroutine GDNcAttrGetDouble(var, name, value, stat, default)
    use gtdata_netcdf_types,   only: GD_NC_VARIABLE, GD_NC_VARIABLE_ENTRY
    use gtdata_netcdf_internal,only: vtable_lookup
    use netcdf_f77, only: nf_noerr, nf_einval, nf_global, nf_char, nf_enomem, nf_inq_att, nf_get_att_Double
    use dc_url,     only: gt_plus
    use gtdata_netcdf_generic, only: get_attr
    use dc_types,   only: STRING, DP
    use dc_string
    implicit none
    type(GD_NC_VARIABLE),       intent(in) :: var
    character(len = *),      intent(in) :: name
    real(DP), intent(out):: value(:)
    integer,            intent(out):: stat
    real(DP),  intent(in), optional:: default

    real(DP),  allocatable:: rbuffer(:)
    character(STRING)      :: cbuffer
    character(STRING), pointer :: lbuffer(:)
    integer           :: attrlen, xtype, i, xferend, iname, varid
    type(GD_NC_VARIABLE_ENTRY):: ent
continue
    stat = vtable_lookup(var, ent)
    if (stat /= nf_noerr) then
        if (present(default)) value(:) = default
        return
    endif
    ! 型と長さを取得
    if (name(1:1) == gt_plus) then
        iname = 2
        varid = nf_global
    else
        iname = 1
        varid = ent%varid
    endif
    stat = nf_inq_att(ent%fileid, varid, name(iname:), xtype=xtype, len=attrlen)
    if (stat /= nf_noerr) then
        if (present(default)) value(:) = default
        return
    endif
    ! 文字型の場合は長さをコンマで分解した語数と読み替える
    if (xtype == nf_char) then
        call get_attr(var, name, cbuffer, "", stat)
        if (stat /= 0) return
        call split(cbuffer, lbuffer, ", ")
        attrlen = size(lbuffer)
    endif
    ! 結果を入れるところがなければ長さだけを伝え終了
    if (size(value) == 0) then
        if (xtype == nf_char) deallocate(lbuffer)
        stat = attrlen
        return
    endif
    ! 型に応じて要求されただけ値を転送
    xferend = min(size(value), attrlen)
    if (present(default)) value(xferend+1: ) = default
    if (xtype == nf_char) then
        do, i = 1, xferend
          value(i) = stod(lbuffer(i))
        enddo
        deallocate(lbuffer)
        stat = attrlen
        return
    else
        allocate(rbuffer(attrlen), stat=stat)
        if (stat /= 0) then
            stat = nf_enomem
            return
        endif
        stat = nf_get_att_Double(ent%fileid, varid, name(iname:), rbuffer)
        if (stat == nf_noerr) then
            value(1:xferend) = rbuffer(1:xferend)
            stat = attrlen
        endif
        deallocate(rbuffer)
        return
    endif
end subroutine GDNcAttrGetDouble
Subroutine :
var :type(GD_NC_VARIABLE), intent(in)
name :character(len = *), intent(in)
value(:) :integer, intent(out)
stat :integer, intent(out)
default :integer, intent(in), optional

[Source]

subroutine GDNcAttrGetInt(var, name, value, stat, default)
    use gtdata_netcdf_types,   only: GD_NC_VARIABLE, GD_NC_VARIABLE_ENTRY
    use gtdata_netcdf_internal,only: vtable_lookup
    use netcdf_f77, only: nf_noerr, nf_einval, nf_global, nf_char, nf_enomem, nf_inq_att, nf_get_att_Int
    use dc_url,     only: gt_plus
    use gtdata_netcdf_generic, only: get_attr
    use dc_types,   only: STRING, DP
    use dc_string
    implicit none
    type(GD_NC_VARIABLE),       intent(in) :: var
    character(len = *),      intent(in) :: name
    integer, intent(out):: value(:)
    integer,            intent(out):: stat
    integer,  intent(in), optional:: default

    integer,  allocatable:: rbuffer(:)
    character(STRING)      :: cbuffer
    character(STRING), pointer :: lbuffer(:)
    integer           :: attrlen, xtype, i, xferend, iname, varid
    type(GD_NC_VARIABLE_ENTRY):: ent
continue
    stat = vtable_lookup(var, ent)
    if (stat /= nf_noerr) then
        if (present(default)) value(:) = default
        return
    endif
    ! 型と長さを取得
    if (name(1:1) == gt_plus) then
        iname = 2
        varid = nf_global
    else
        iname = 1
        varid = ent%varid
    endif
    stat = nf_inq_att(ent%fileid, varid, name(iname:), xtype=xtype, len=attrlen)
    if (stat /= nf_noerr) then
        if (present(default)) value(:) = default
        return
    endif
    ! 文字型の場合は長さをコンマで分解した語数と読み替える
    if (xtype == nf_char) then
        call get_attr(var, name, cbuffer, "", stat)
        if (stat /= 0) return
        call split(cbuffer, lbuffer, ", ")
        attrlen = size(lbuffer)
    endif
    ! 結果を入れるところがなければ長さだけを伝え終了
    if (size(value) == 0) then
        if (xtype == nf_char) deallocate(lbuffer)
        stat = attrlen
        return
    endif
    ! 型に応じて要求されただけ値を転送
    xferend = min(size(value), attrlen)
    if (present(default)) value(xferend+1: ) = default
    if (xtype == nf_char) then
        do, i = 1, xferend
          value(i) = stod(lbuffer(i))
        enddo
        deallocate(lbuffer)
        stat = attrlen
        return
    else
        allocate(rbuffer(attrlen), stat=stat)
        if (stat /= 0) then
            stat = nf_enomem
            return
        endif
        stat = nf_get_att_Int(ent%fileid, varid, name(iname:), rbuffer)
        if (stat == nf_noerr) then
            value(1:xferend) = rbuffer(1:xferend)
            stat = attrlen
        endif
        deallocate(rbuffer)
        return
    endif
end subroutine GDNcAttrGetInt
Subroutine :
var :type(GD_NC_VARIABLE), intent(in)
name :character(len = *), intent(in)
value(:) :real, intent(out)
stat :integer, intent(out)
default :real, intent(in), optional

[Source]

subroutine GDNcAttrGetReal(var, name, value, stat, default)
    use gtdata_netcdf_types,   only: GD_NC_VARIABLE, GD_NC_VARIABLE_ENTRY
    use gtdata_netcdf_internal,only: vtable_lookup
    use netcdf_f77, only: nf_noerr, nf_einval, nf_global, nf_char, nf_enomem, nf_inq_att, nf_get_att_Real
    use dc_url,     only: gt_plus
    use gtdata_netcdf_generic, only: get_attr
    use dc_types,   only: STRING, DP
    use dc_string
    implicit none
    type(GD_NC_VARIABLE),       intent(in) :: var
    character(len = *),      intent(in) :: name
    real, intent(out):: value(:)
    integer,            intent(out):: stat
    real,  intent(in), optional:: default

    real,  allocatable:: rbuffer(:)
    character(STRING)      :: cbuffer
    character(STRING), pointer :: lbuffer(:)
    integer           :: attrlen, xtype, i, xferend, iname, varid
    type(GD_NC_VARIABLE_ENTRY):: ent
continue
    stat = vtable_lookup(var, ent)
    if (stat /= nf_noerr) then
        if (present(default)) value(:) = default
        return
    endif
    ! 型と長さを取得
    if (name(1:1) == gt_plus) then
        iname = 2
        varid = nf_global
    else
        iname = 1
        varid = ent%varid
    endif
    stat = nf_inq_att(ent%fileid, varid, name(iname:), xtype=xtype, len=attrlen)
    if (stat /= nf_noerr) then
        if (present(default)) value(:) = default
        return
    endif
    ! 文字型の場合は長さをコンマで分解した語数と読み替える
    if (xtype == nf_char) then
        call get_attr(var, name, cbuffer, "", stat)
        if (stat /= 0) return
        call split(cbuffer, lbuffer, ", ")
        attrlen = size(lbuffer)
    endif
    ! 結果を入れるところがなければ長さだけを伝え終了
    if (size(value) == 0) then
        if (xtype == nf_char) deallocate(lbuffer)
        stat = attrlen
        return
    endif
    ! 型に応じて要求されただけ値を転送
    xferend = min(size(value), attrlen)
    if (present(default)) value(xferend+1: ) = default
    if (xtype == nf_char) then
        do, i = 1, xferend
          value(i) = stod(lbuffer(i))
        enddo
        deallocate(lbuffer)
        stat = attrlen
        return
    else
        allocate(rbuffer(attrlen), stat=stat)
        if (stat /= 0) then
            stat = nf_enomem
            return
        endif
        stat = nf_get_att_Real(ent%fileid, varid, name(iname:), rbuffer)
        if (stat == nf_noerr) then
            value(1:xferend) = rbuffer(1:xferend)
            stat = attrlen
        endif
        deallocate(rbuffer)
        return
    endif
end subroutine GDNcAttrGetReal