gdncvarcreate.f90 File Reference

Go to the source code of this file.

Functions/Subroutines

subroutine gdncvarcreate (var, url, xtype, dims, overwrite, err)
 

Function/Subroutine Documentation

◆ gdncvarcreate()

subroutine gdncvarcreate ( type(gd_nc_variable), intent(out)  var,
character(len = *), intent(in)  url,
character(len = *), intent(in)  xtype,
type(gd_nc_variable), dimension(:), intent(in)  dims,
logical, intent(in), optional  overwrite,
logical, intent(out), optional  err 
)

Definition at line 14 of file gdncvarcreate.f90.

References dc_trace::beginsub(), dc_trace::dbgmessage(), dc_trace::endsub(), gdncfiledefinemode(), gdncfileopen(), dc_error::gt_edimmultidim, dc_error::gt_edimnodim, dc_error::gt_enomem, dc_error::gt_eotherfile, dc_error::storeerror(), dc_types::string, gtdata_netcdf_internal::vtable_add(), and gtdata_netcdf_internal::vtable_lookup().

14  !
15  !== 変数作成
16  !
17  ! 変数 URL *url* に変数を作成します.
18  ! 変数が依存する次元を *dims* に与えます.
19  ! 返される引数 *var* には変数 ID などの情報が格納されます.
20  !
21  ! *overwrite* に .true. を設定すると上書き可能なモードになります.
22  ! デフォルトは上書き不可です.
23  ! *err* を与える場合, 次元変数生成時にエラーが生じても
24  ! プログラムを終了せず, *err* に .false. が返ります.
25  !
27  use dc_types, only: string
28  use dc_string, only: strieq
31  use dc_url, only: urlsplit
32  use dc_trace, only: beginsub, endsub, dbgmessage
33  use gtdata_netcdf_generic, only: tostring ! for debug
34  use netcdf, only: &
35  & nf90_noerr, nf90_float, nf90_double, nf90_int, nf90_char, nf90_ebaddim, nf90_def_var
38  implicit none
39  type(gd_nc_variable), intent(out):: var
40  character(len = *), intent(in):: url
41  character(len = *), intent(in):: xtype
42  type(gd_nc_variable), intent(in):: dims(:)
43  logical, intent(in), optional:: overwrite
44  logical, intent(out), optional:: err
45  type(gd_nc_variable_search):: ent
46  type(gd_nc_variable_entry):: ent_dim
47  character(len = string):: filename, varname
48  integer, allocatable:: dimids(:)
49  integer:: stat, nvdims, i
50  integer:: nc_xtype
51  logical:: clobber
52  intrinsic trim
53  character(len = *), parameter:: subnam = "GDNcVarCreate"
54 continue
55  clobber = .false.
56  if (present(overwrite)) clobber = overwrite
57  call beginsub(subnam)
58  call dbgmessage('url=%c', c1=trim(url))
59  call dbgmessage('xtype=%c', c1=trim(xtype))
60  call dbgmessage('dims=(/%*d/)', i=(/dims(:)%id/), n=(/size(dims)/))
61  call dbgmessage('ovwr=%y', l=(/clobber/))
62 
63  ! もし必要ならファイル作成
64  call urlsplit(url, filename, varname)
65  call gdncfileopen(ent%fileid, filename, stat=stat, writable=.true., &
66  & overwrite=clobber)
67  if (stat /= nf90_noerr) goto 999
68 
69  ! 次元にまつわる準備
70  nvdims = size(dims)
71  allocate(dimids(max(1, nvdims)), stat=stat)
72  if (stat /= 0) then
73  stat = gt_enomem
74  goto 999
75  end if
76  do, i = 1, nvdims
77  stat = vtable_lookup(dims(i), ent_dim)
78  if (stat /= nf90_noerr) then
79  stat = nf90_ebaddim
80  goto 999
81  endif
82  if (ent%fileid /= ent_dim%fileid) then
83  stat = gt_eotherfile
84  goto 999
85  endif
86  if (ent_dim%dimid <= 0) then
87  stat = gt_edimmultidim
88  goto 999
89  endif
90  dimids(i) = ent_dim%dimid
91  enddo
92  ent%dimid = 0
93 
94  ! 変数の型の判定
95  nc_xtype = nf90_float
96  if (strieq(xtype, "double") .or. strieq(xtype, "DOUBLEPRECISION")) then
97  nc_xtype = nf90_double
98  endif
99  if (strieq(xtype, "int") .or. strieq(xtype, "INTEGER")) then
100  nc_xtype = nf90_int
101  endif
102  if (strieq(xtype, "char") .or. strieq(xtype, "CHARACTER")) then
103  nc_xtype = nf90_char
104  endif
105 
106  ! 本当の変数作成操作
107  stat = gdncfiledefinemode(ent%fileid)
108  if (stat /= nf90_noerr) goto 999
109  if ( nvdims == 0 ) then
110  stat = nf90_def_var(ent%fileid, name = trim(varname), &
111  & xtype = nc_xtype, varid=ent%varid)
112  else
113  stat = nf90_def_var(ent%fileid, name = trim(varname), &
114  & xtype = nc_xtype, dimids = dimids, varid=ent%varid)
115  end if
116  if (stat /= nf90_noerr) goto 999
117 
118  ! 登録
119  stat = vtable_add(var, ent)
120 
121 999 continue
122  if (allocated(dimids)) deallocate(dimids)
123  if (stat /= nf90_noerr) var % id = -1
124  call storeerror(stat, subnam, err, cause_c=url)
125  call endsub(subnam, 'stat=%d, var.id=%d', i=(/stat, var % id/))
integer function, public vtable_lookup(var, entry)
integer function gdncfiledefinemode(fileid)
integer, parameter, public gt_edimnodim
Definition: dc_error.f90:529
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, parameter, public gt_eotherfile
Definition: dc_error.f90:535
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
integer, parameter, public gt_edimmultidim
Definition: dc_error.f90:530
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
Here is the call graph for this function: