anattrgetnum.f90

Path: src/anattrgetnum.f90
Last Update: Thu Sep 08 22:21:48 JST 2005
    Copyright (C) GFD Dennou Club, 2000-2005. All rights reserved.

Get AN_VARIABLES

This file is created by "anattergettype.m4" by m4 command using "intrinsic_types.m4". Don‘t edit each files directly.

お客様向きではないけれど、情報落ちのないインターフェイスということで.… stat < 0: エラー、あるいはその属性は存在しなかった stat = 0 … size(value): その属性を全部読み取った。サイズは stat 個 stat > size(value): 配列長不足のため属性が全部読み取れなかった。

                     サイズは stat 個必要

バグ:

 属性が文字型で STRING 文字を越える場合、GT_ECHARSHORT が返る

Methods

Included Modules

an_types an_vartable netcdf_f77 dc_url an_generic dc_types dc_string

Public Instance methods

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

[Source]



subroutine ANAttrGetDouble(var, name, value, stat, default)

    implicit none
    type(AN_VARIABLE),       intent(in) :: var
    character(len = *),      intent(in) :: name
    real(DP), intent(out):: value(:)
    integer(INTK),            intent(out):: stat
    real(DP),  intent(in), optional:: default

    real(DP),  allocatable:: rbuffer(:)
    character(STRING)      :: cbuffer
    type(STRING_LIST)      :: lbuffer
    integer(INTK)           :: attrlen, xtype, i, xferend, iname, varid
    type(AN_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(lbuffer, cbuffer, ", ")
        attrlen = len(lbuffer)
    endif
    ! 結果を入れるところがなければ長さだけを伝え終了
    if (size(value) == 0) then
        if (xtype == nf_char) call dispose(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(element(lbuffer, i))
        enddo
        call dispose(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 ANAttrGetDouble
var :type(AN_VARIABLE), intent(in)
name :character(len = *), intent(in)
value(:) :integer(INTK), intent(out)
stat :integer(INTK), intent(out)
default :integer(INTK), intent(in), optional

[Source]


subroutine ANAttrGetInt(var, name, value, stat, default)

    implicit none
    type(AN_VARIABLE),       intent(in) :: var
    character(len = *),      intent(in) :: name
    integer(INTK), intent(out):: value(:)
    integer(INTK),            intent(out):: stat
    integer(INTK),  intent(in), optional:: default

    integer(INTK),  allocatable:: rbuffer(:)
    character(STRING)      :: cbuffer
    type(STRING_LIST)      :: lbuffer
    integer(INTK)           :: attrlen, xtype, i, xferend, iname, varid
    type(AN_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(lbuffer, cbuffer, ", ")
        attrlen = len(lbuffer)
    endif
    ! 結果を入れるところがなければ長さだけを伝え終了
    if (size(value) == 0) then
        if (xtype == nf_char) call dispose(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(element(lbuffer, i))
        enddo
        call dispose(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 ANAttrGetInt
var :type(AN_VARIABLE), intent(in)
name :character(len = *), intent(in)
value(:) :real(SP), intent(out)
stat :integer(INTK), intent(out)
default :real(SP), intent(in), optional

[Source]



subroutine ANAttrGetReal(var, name, value, stat, default)

    implicit none
    type(AN_VARIABLE),       intent(in) :: var
    character(len = *),      intent(in) :: name
    real(SP), intent(out):: value(:)
    integer(INTK),            intent(out):: stat
    real(SP),  intent(in), optional:: default

    real(SP),  allocatable:: rbuffer(:)
    character(STRING)      :: cbuffer
    type(STRING_LIST)      :: lbuffer
    integer(INTK)           :: attrlen, xtype, i, xferend, iname, varid
    type(AN_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(lbuffer, cbuffer, ", ")
        attrlen = len(lbuffer)
    endif
    ! 結果を入れるところがなければ長さだけを伝え終了
    if (size(value) == 0) then
        if (xtype == nf_char) call dispose(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(element(lbuffer, i))
        enddo
        call dispose(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 ANAttrGetReal

[Validate]