gdncvarcreated.f90

Path: gtdata/gtdata_netcdf/gdncvarcreated.f90
Last Update: Mon May 25 18:51:59 +0900 2009

netCDF ファイルへ次元変数作成

Authors:Eizi TOYODA, Yasuhiro MORIKAWA
Version:$Id: gdncvarcreated.f90,v 1.2 2009-05-25 09:51:59 morikawa Exp $
Tag Name:$Name: gtool5-20090809 $
Copyright:Copyright (C) GFD Dennou Club, 2000-2005. All rights reserved.
License:See COPYRIGHT

以下のサブルーチン、関数は gtdata_netcdf_generic から gtdata_netcdf_generic#Create として提供されます。

Methods

Included Modules

gtdata_netcdf_types gtdata_netcdf_internal dc_string dc_types dc_url dc_trace netcdf_f77 gtdata_netcdf_file_generic dc_error

Public Instance methods

Subroutine :
var :type(GD_NC_VARIABLE), intent(out)
url :character(len = *), intent(in)
xtype :character(len = *), intent(in)
length :integer, intent(in)
overwrite :logical, intent(in), optional
err :logical, intent(out), optional

次元変数作成

変数 URL url に次元変数を作成します. 次元変数の長さを length に与えます. 返される引数 var には変数 ID などの情報が格納されます.

overwrite に .true. を設定すると上書き可能なモードになります. デフォルトは上書き不可です. err を与える場合, 次元変数生成時にエラーが生じても プログラムを終了せず, err に .false. が返ります.

[Source]

subroutine GDNcVarCreateD(var, url, xtype, length, overwrite, err)
  !
  !== 次元変数作成
  !
  ! 変数 URL *url* に次元変数を作成します.
  ! 次元変数の長さを *length* に与えます.
  ! 返される引数 *var* には変数 ID などの情報が格納されます.
  !
  ! *overwrite* に .true. を設定すると上書き可能なモードになります.
  ! デフォルトは上書き不可です.
  ! *err* を与える場合, 次元変数生成時にエラーが生じても
  ! プログラムを終了せず, *err* に .false. が返ります.
  !
  use gtdata_netcdf_types, only: GD_NC_VARIABLE, GD_NC_VARIABLE_SEARCH
  use gtdata_netcdf_internal, only: vtable_add
  use dc_string, only: strieq
  use dc_types, only: string
  use dc_url, only: UrlSplit
  use dc_trace, only: BeginSub, EndSub, DbgMessage
  use netcdf_f77, only: nf_noerr, nf_real, nf_int, nf_double, nf_def_var, nf_def_dim
  use gtdata_netcdf_file_generic, only: GDNcFileOpen, GDNcFileDefineMode
  use dc_error, only: StoreError, gt_enomem
  implicit none
  type(GD_NC_VARIABLE), intent(out):: var
  character(len = *), intent(in):: url
  character(len = *), intent(in):: xtype
  integer, intent(in):: length
  logical, intent(in), optional:: overwrite
  logical, intent(out), optional:: err
  type(GD_NC_VARIABLE_SEARCH):: ent
  character(len = string):: filename, varname, cause_c
  integer:: stat
  integer:: nc_xtype
  character(len = *), parameter:: subname = "GDNcVarCreateD"
continue
  call BeginSub(subname, 'url=<%c>, xtype=<%c>, length=<%d>', c1=trim(url), c2=trim(xtype), i=(/length/))
  cause_c = trim(url)
  !
  ! --- ファイルを用意 ---
  call UrlSplit(url, file=filename, var=varname)
  call GDNcFileOpen(ent%fileid, filename, stat=stat, writable=.TRUE., overwrite=overwrite)
  if (stat /= NF_NOERR) goto 999
  stat = GDNcFileDefineMode(ent%fileid)
  if (stat /= NF_NOERR) goto 999
  !
  ! --- 型の決定 ---
  nc_xtype = NF_REAL
  if (strieq(xtype, "double") .or. strieq(xtype, "DOUBLEPRECISION")) then
    nc_xtype = NF_DOUBLE
  endif
  if (strieq(xtype, "int") .or. strieq(xtype, "INTEGER")) then
    nc_xtype = NF_INT
  endif
  !
  ! --- 次元変数の作成 ---
  stat = nf_def_dim(ent%fileid, trim(varname), len=length, dimid=ent%dimid)
  if (stat /= NF_NOERR) goto 999
  stat = nf_def_var(ent%fileid, trim(varname), xtype=nc_xtype, ndims=1, dimids=(/ent%dimid/), varid=ent%varid)
  if (stat /= NF_NOERR) goto 999
  !
  stat = vtable_add(var, ent)
  if (stat /= NF_NOERR) goto 999

999 continue
  call StoreError(stat, subname, err, cause_c=cause_c)
  if (stat /= NF_NOERR) var = GD_NC_VARIABLE(-1)
  call EndSub(subname, 'stat=%d', i=(/stat/))
end subroutine GDNcVarCreateD

[Validate]