! -*- coding: utf-8; mode: f90 -*-
!-------------------------------------------------------------------------------------
! Copyright (c) 2000-2016 Gtool Development Group. All rights reserved.
!-------------------------------------------------------------------------------------
! ** Important**
!
! This file is generated from gdncattrgetnum.erb by ERB included Ruby 2.7.4.
! Please do not edit this file directly. @see "gdncattrgetnum.erb"
!-------------------------------------------------------------------------------------
!
!
!= Get GD_NC_VARIABLES
!
! Authors::   Eizi TOYODA, Yasuhiro MORIKAWA
! Version::   $Id: gdncattrgetnum.rb2f90,v 1.2 2009-05-25 09:51:59 morikawa Exp $
! Tag Name::  $Name:  $
! Copyright:: Copyright (C) GFD Dennou Club, 2000-2009. All rights reserved.
! License::   See COPYRIGHT[link:../../COPYRIGHT]
!
! お客様向きではないけれど、情報落ちのないインターフェイスということで....
!
! stat < 0                :: エラー、あるいはその属性は存在しなかった
! stat = 0 ... size(value):: その属性を全部読み取った。サイズは stat 個
! stat > size(value)      :: 配列長不足のため属性が全部読み取れなかった。
!                            サイズは stat 個必要
!
! * バグ:
!   * 属性が文字型で STRING 文字を越える場合、GT_ECHARSHORT が返る
!
!
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, only: &
    & NF90_NOERR, &
    & NF90_GLOBAL, &
    & NF90_CHAR, &
    & NF90_ENOMEM, &
    & NF90_INQUIRE_ATTRIBUTE, &
    & NF90_GET_ATT
  use dc_url,     only: gt_plus
  use gtdata_netcdf_generic, only: get_attr
  use dc_types,   only: STRING
  use dc_string
  implicit none
  type(GD_NC_VARIABLE),       intent(in) :: var
  character(len = *),      intent(in) :: name
  integer,            intent(out):: stat
  integer, intent(out):: value(:)
  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 /= NF90_NOERR) then
    if (present(default)) value(:) = default
    return
  endif
  ! 型と長さを取得
  if (name(1:1) == gt_plus) then
    iname = 2
    varid = NF90_GLOBAL
  else
    iname = 1
    varid = ent%varid
  endif
  stat = NF90_INQUIRE_ATTRIBUTE(ent%fileid, varid, &
    name = name(iname:), xtype = xtype, len = attrlen)
  if (stat /= NF90_NOERR) then
    if (present(default)) value(:) = default
    return
  endif
  ! 文字型の場合は長さをコンマで分解した語数と読み替える
  if (xtype == NF90_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 == NF90_CHAR) deallocate(lbuffer)
    stat = attrlen
    return
  endif
  ! 型に応じて要求されただけ値を転送
  xferend = min(size(value), attrlen)
  if (present(default)) value(xferend+1: ) = default
  if (xtype == NF90_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 = NF90_ENOMEM
      return
    endif
    stat = NF90_GET_ATT(ent%fileid, varid, name(iname:), rbuffer)
    if (stat == NF90_NOERR) then
      value(1:xferend) = rbuffer(1:xferend)
      stat = attrlen
    endif
    deallocate(rbuffer)
    return
  endif
end subroutine GDNcAttrGetInt
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, only: &
    & NF90_NOERR, &
    & NF90_GLOBAL, &
    & NF90_CHAR, &
    & NF90_ENOMEM, &
    & NF90_INQUIRE_ATTRIBUTE, &
    & NF90_GET_ATT
  use dc_url,     only: gt_plus
  use gtdata_netcdf_generic, only: get_attr
  use dc_types,   only: STRING, SP
  use dc_string
  implicit none
  type(GD_NC_VARIABLE),       intent(in) :: var
  character(len = *),      intent(in) :: name
  integer,            intent(out):: stat
  real(SP), intent(out):: value(:)
  real(SP),  intent(in), optional:: default
  real(SP),  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 /= NF90_NOERR) then
    if (present(default)) value(:) = default
    return
  endif
  ! 型と長さを取得
  if (name(1:1) == gt_plus) then
    iname = 2
    varid = NF90_GLOBAL
  else
    iname = 1
    varid = ent%varid
  endif
  stat = NF90_INQUIRE_ATTRIBUTE(ent%fileid, varid, &
    name = name(iname:), xtype = xtype, len = attrlen)
  if (stat /= NF90_NOERR) then
    if (present(default)) value(:) = default
    return
  endif
  ! 文字型の場合は長さをコンマで分解した語数と読み替える
  if (xtype == NF90_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 == NF90_CHAR) deallocate(lbuffer)
    stat = attrlen
    return
  endif
  ! 型に応じて要求されただけ値を転送
  xferend = min(size(value), attrlen)
  if (present(default)) value(xferend+1: ) = default
  if (xtype == NF90_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 = NF90_ENOMEM
      return
    endif
    stat = NF90_GET_ATT(ent%fileid, varid, name(iname:), rbuffer)
    if (stat == NF90_NOERR) then
      value(1:xferend) = rbuffer(1:xferend)
      stat = attrlen
    endif
    deallocate(rbuffer)
    return
  endif
end subroutine GDNcAttrGetReal
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, only: &
    & NF90_NOERR, &
    & NF90_GLOBAL, &
    & NF90_CHAR, &
    & NF90_ENOMEM, &
    & NF90_INQUIRE_ATTRIBUTE, &
    & NF90_GET_ATT
  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):: stat
  real(DP), intent(out):: value(:)
  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 /= NF90_NOERR) then
    if (present(default)) value(:) = default
    return
  endif
  ! 型と長さを取得
  if (name(1:1) == gt_plus) then
    iname = 2
    varid = NF90_GLOBAL
  else
    iname = 1
    varid = ent%varid
  endif
  stat = NF90_INQUIRE_ATTRIBUTE(ent%fileid, varid, &
    name = name(iname:), xtype = xtype, len = attrlen)
  if (stat /= NF90_NOERR) then
    if (present(default)) value(:) = default
    return
  endif
  ! 文字型の場合は長さをコンマで分解した語数と読み替える
  if (xtype == NF90_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 == NF90_CHAR) deallocate(lbuffer)
    stat = attrlen
    return
  endif
  ! 型に応じて要求されただけ値を転送
  xferend = min(size(value), attrlen)
  if (present(default)) value(xferend+1: ) = default
  if (xtype == NF90_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 = NF90_ENOMEM
      return
    endif
    stat = NF90_GET_ATT(ent%fileid, varid, name(iname:), rbuffer)
    if (stat == NF90_NOERR) then
      value(1:xferend) = rbuffer(1:xferend)
      stat = attrlen
    endif
    deallocate(rbuffer)
    return
  endif
end subroutine GDNcAttrGetDouble
