anvaropenbydimord.f90

Path: anvaropenbydimord.f90
Last Update: Wed Jul 20 18:22:22 JST 2005

anvaropenbydimord.f90 - open(AN_VARIABLE, AN_VARIABLE, integer) Copyright (C) GFD Dennou Club, 2000. All rights reserved

open(dimvar, var, dimord, [err]) は 既に開かれた変数 var の ord 番目の次元にあたる変数を 開き dimvar に格納する。

Methods

Included Modules

an_types an_vartable an_file netcdf_f77 dc_error

Public Instance methods

Subroutine :
var :type(AN_VARIABLE), intent(out)
src_var :type(AN_VARIABLE), intent(in)
dimord :integer, intent(in)
err :logical, intent(out), optional

[Source]

subroutine ANVarOpenByDimOrd(var, src_var, dimord, err)
    use an_types, only: AN_VARIABLE, an_variable_entry, an_variable_search
    use an_vartable, only: vtable_lookup, vtable_add
    use an_file
    use netcdf_f77
    use dc_error
    implicit none
    type(AN_VARIABLE), intent(out):: var
    type(AN_VARIABLE), intent(in):: src_var
    integer, intent(in):: dimord
    logical, intent(out), optional:: err
    type(an_variable_entry):: src_ent
    type(an_variable_search):: ent
    character(len = NF_MAX_NAME):: dimname
    integer:: stat
continue
    stat = vtable_lookup(src_var, src_ent)
    if (stat /= NF_NOERR) goto 999
    !
    if (dimord <= 0) then
        !
        ! dimord == 0 の特別扱いは廃止 (gtdata 層でやる)。
        ! (同じハンドルを返すと参照カウントが増えないので多重クローズが
        ! 禁止されてしまい、Open のセマンティクスに不適合なため。)
        !
        var = src_var
        stat = NF_EINVAL
        goto 999
    endif
    if (.not. associated(src_ent%dimids)) then
        stat = GT_ENOMOREDIMS
        goto 999
    else if (dimord > size(src_ent%dimids)) then
        stat = GT_ENOMOREDIMS
        goto 999
    endif
    !
    ! 決定された次元 id を使って変数 id の確定を試みる
    !
    ent%fileid = src_ent%fileid
    ent%dimid = src_ent%dimids(dimord)
    stat = nf_inq_dimname(ent%fileid, ent%dimid, dimname)
    if (stat /= NF_NOERR) goto 999
    stat = nf_inq_varid(ent%fileid, dimname, ent%varid)
    if (stat == NF_ENOTVAR) then
        ! 変数が対応しなければ無視するだけ
        ent%varid = 0
        stat = NF_NOERR
    else if (stat /= NF_NOERR) then
        goto 999
    endif
    !
    ! 成功しそうなので参照カウント値を増やす
    call ANFileReopen(ent%fileid)
    ! 不適合な変数が開かれている場合は vtable_add が変数部を無視する。
    stat = vtable_add(var, ent)
    if (stat /= nf_noerr) then
        call anfileclose(ent%fileid)
        goto 999
    endif
    call StoreError(nf_noerr, 'ANVarOpenByDimOrd', err)
    return

999 continue
    call StoreError(stat, 'ANVarOpenByDimOrd', err, cause_i=dimord)
    return
end subroutine

[Validate]