gdncvarcreated.f90
Go to the documentation of this file.
1 !
2 != netCDF ファイルへ次元変数作成
3 !
4 ! Authors:: Eizi TOYODA, Yasuhiro MORIKAWA
5 ! Version:: $Id: gdncvarcreated.f90,v 1.3 2009-10-11 07:36:34 morikawa Exp $
6 ! Tag Name:: $Name: $
7 ! Copyright:: Copyright (C) GFD Dennou Club, 2000-2005. All rights reserved.
8 ! License:: See COPYRIGHT[link:../../COPYRIGHT]
9 !
10 ! 以下のサブルーチン、関数は gtdata_netcdf_generic から gtdata_netcdf_generic#Create
11 ! として提供されます。
12 
13 subroutine gdncvarcreated(var, url, xtype, length, overwrite, err)
14  !
15  !== 次元変数作成
16  !
17  ! 変数 URL *url* に次元変数を作成します.
18  ! 次元変数の長さを *length* に与えます.
19  ! 返される引数 *var* には変数 ID などの情報が格納されます.
20  !
21  ! *overwrite* に .true. を設定すると上書き可能なモードになります.
22  ! デフォルトは上書き不可です.
23  ! *err* を与える場合, 次元変数生成時にエラーが生じても
24  ! プログラムを終了せず, *err* に .false. が返ります.
25  !
28  use dc_string, only: strieq
29  use dc_types, only: string
30  use dc_url, only: urlsplit
31  use dc_trace, only: beginsub, endsub, dbgmessage
32  use netcdf, only: nf90_noerr, nf90_float, nf90_int, nf90_double, nf90_char, &
33  & nf90_def_var, nf90_def_dim
35  use dc_error, only: storeerror, gt_enomem
36  implicit none
37  type(gd_nc_variable), intent(out):: var
38  character(len = *), intent(in):: url
39  character(len = *), intent(in):: xtype
40  integer, intent(in):: length
41  logical, intent(in), optional:: overwrite
42  logical, intent(out), optional:: err
43  type(gd_nc_variable_search):: ent
44  character(len = string):: filename, varname, cause_c
45  integer:: stat
46  integer:: nc_xtype
47  character(len = *), parameter:: subname = "GDNcVarCreateD"
48 continue
49  call beginsub(subname, 'url=<%c>, xtype=<%c>, length=<%d>', &
50  & c1=trim(url), c2=trim(xtype), i=(/length/))
51  cause_c = trim(url)
52  !
53  ! --- ファイルを用意 ---
54  call urlsplit(url, file=filename, var=varname)
55  call gdncfileopen(ent%fileid, filename, stat=stat, writable=.true., &
56  & overwrite=overwrite)
57  if (stat /= nf90_noerr) goto 999
58  stat = gdncfiledefinemode(ent%fileid)
59  if (stat /= nf90_noerr) goto 999
60  !
61  ! --- 型の決定 ---
62  nc_xtype = nf90_float
63  if (strieq(xtype, "double") .or. strieq(xtype, "DOUBLEPRECISION")) then
64  nc_xtype = nf90_double
65  endif
66  if (strieq(xtype, "int") .or. strieq(xtype, "INTEGER")) then
67  nc_xtype = nf90_int
68  endif
69  if (strieq(xtype, "char") .or. strieq(xtype, "CHARACTER")) then
70  nc_xtype = nf90_char
71  endif
72  !
73  ! --- 次元変数の作成 ---
74  stat = nf90_def_dim(ent%fileid, trim(varname), len=length, dimid=ent%dimid)
75  if (stat /= nf90_noerr) goto 999
76  stat = nf90_def_var(ent%fileid, trim(varname), &
77  & xtype=nc_xtype, dimids=(/ent%dimid/), varid=ent%varid)
78  if (stat /= nf90_noerr) goto 999
79  !
80  stat = vtable_add(var, ent)
81  if (stat /= nf90_noerr) goto 999
82 
83 999 continue
84  call storeerror(stat, subname, err, cause_c=cause_c)
85  if (stat /= nf90_noerr) var = gd_nc_variable(-1)
86  call endsub(subname, 'stat=%d', i=(/stat/))
87 end subroutine gdncvarcreated
integer function gdncfiledefinemode(fileid)
subroutine gdncvarcreated(var, url, xtype, length, overwrite, err)
integer, parameter, public gt_enomem
Definition: dc_error.f90:534
subroutine, public storeerror(number, where, err, cause_c, cause_i)
Definition: dc_error.f90:830
integer function, public vtable_add(var, entry)
subroutine, public dbgmessage(fmt, i, r, d, L, n, c1, c2, c3, ca)
Definition: dc_trace.f90:509
subroutine, public beginsub(name, fmt, i, r, d, L, n, c1, c2, c3, ca, version)
Definition: dc_trace.f90:351
subroutine gdncfileopen(fileid, filename, writable, overwrite, stat, err)
Definition: gdncfileopen.f90:2
文字型変数の操作.
Definition: dc_string.f90:24
種別型パラメタを提供します。
Definition: dc_types.f90:49
subroutine, public endsub(name, fmt, i, r, d, L, n, c1, c2, c3, ca)
Definition: dc_trace.f90:446
integer, parameter, public string
文字列を保持する 文字型変数の種別型パラメタ
Definition: dc_types.f90:118