gdncattrgetnum.f90 File Reference

Go to the source code of this file.

Functions/Subroutines

subroutine gdncattrgetint (var, name, value, stat, default)
 
subroutine gdncattrgetreal (var, name, value, stat, default)
 
subroutine gdncattrgetdouble (var, name, value, stat, default)
 

Function/Subroutine Documentation

◆ gdncattrgetdouble()

subroutine gdncattrgetdouble ( type(gd_nc_variable), intent(in)  var,
character(len = *), intent(in)  name,
real(dp), dimension(:), intent(out)  value,
integer, intent(out)  stat,
real(dp), intent(in), optional  default 
)

Definition at line 200 of file gdncattrgetnum.f90.

References dc_types::dp, dc_url::gt_plus, dc_types::string, and gtdata_netcdf_internal::vtable_lookup().

202  use netcdf, only: &
203  & nf90_noerr, &
204  & nf90_global, &
205  & nf90_char, &
206  & nf90_enomem, &
207  & nf90_inquire_attribute, &
208  & nf90_get_att
209  use dc_url, only: gt_plus
211  use dc_types, only: string, dp
212  use dc_string
213  implicit none
214  type(gd_nc_variable), intent(in) :: var
215  character(len = *), intent(in) :: name
216  integer, intent(out):: stat
217  real(DP), intent(out):: value(:)
218  real(DP), intent(in), optional:: default
219  real(DP), allocatable:: rbuffer(:)
220  character(STRING) :: cbuffer
221  character(STRING), pointer :: lbuffer(:)
222  integer :: attrlen, xtype, i, xferend, iname, varid
223  type(gd_nc_variable_entry):: ent
224  continue
225  stat = vtable_lookup(var, ent)
226  if (stat /= nf90_noerr) then
227  if (present(default)) value(:) = default
228  return
229  endif
230  ! 型と長さを取得
231  if (name(1:1) == gt_plus) then
232  iname = 2
233  varid = nf90_global
234  else
235  iname = 1
236  varid = ent%varid
237  endif
238  stat = nf90_inquire_attribute(ent%fileid, varid, &
239  name = name(iname:), xtype = xtype, len = attrlen)
240  if (stat /= nf90_noerr) then
241  if (present(default)) value(:) = default
242  return
243  endif
244  ! 文字型の場合は長さをコンマで分解した語数と読み替える
245  if (xtype == nf90_char) then
246  call get_attr(var, name, cbuffer, "", stat)
247  if (stat /= 0) return
248  call split(cbuffer, lbuffer, ", ")
249  attrlen = size(lbuffer)
250  endif
251  ! 結果を入れるところがなければ長さだけを伝え終了
252  if (size(value) == 0) then
253  if (xtype == nf90_char) deallocate(lbuffer)
254  stat = attrlen
255  return
256  endif
257  ! 型に応じて要求されただけ値を転送
258  xferend = min(size(value), attrlen)
259  if (present(default)) value(xferend+1: ) = default
260  if (xtype == nf90_char) then
261  do, i = 1, xferend
262  value(i) = stod(lbuffer(i))
263  enddo
264  deallocate(lbuffer)
265  stat = attrlen
266  return
267  else
268  allocate(rbuffer(attrlen), stat=stat)
269  if (stat /= 0) then
270  stat = nf90_enomem
271  return
272  endif
273  stat = nf90_get_att(ent%fileid, varid, name(iname:), rbuffer)
274  if (stat == nf90_noerr) then
275  value(1:xferend) = rbuffer(1:xferend)
276  stat = attrlen
277  endif
278  deallocate(rbuffer)
279  return
280  endif
integer function, public vtable_lookup(var, entry)
character, parameter, public gt_plus
Definition: dc_url.f90:92
integer, parameter, public dp
倍精度実数型変数
Definition: dc_types.f90:83
文字型変数の操作.
Definition: dc_string.f90:24
種別型パラメタを提供します。
Definition: dc_types.f90:49
integer, parameter, public string
文字列を保持する 文字型変数の種別型パラメタ
Definition: dc_types.f90:118
Here is the call graph for this function:

◆ gdncattrgetint()

subroutine gdncattrgetint ( type(gd_nc_variable), intent(in)  var,
character(len = *), intent(in)  name,
integer, dimension(:), intent(out)  value,
integer, intent(out)  stat,
integer, intent(in), optional  default 
)

Definition at line 32 of file gdncattrgetnum.f90.

References dc_url::gt_plus, dc_types::string, and gtdata_netcdf_internal::vtable_lookup().

34  use netcdf, only: &
35  & nf90_noerr, &
36  & nf90_global, &
37  & nf90_char, &
38  & nf90_enomem, &
39  & nf90_inquire_attribute, &
40  & nf90_get_att
41  use dc_url, only: gt_plus
43  use dc_types, only: string
44  use dc_string
45  implicit none
46  type(gd_nc_variable), intent(in) :: var
47  character(len = *), intent(in) :: name
48  integer, intent(out):: stat
49  integer, intent(out):: value(:)
50  integer, intent(in), optional:: default
51  integer, allocatable:: rbuffer(:)
52  character(STRING) :: cbuffer
53  character(STRING), pointer :: lbuffer(:)
54  integer :: attrlen, xtype, i, xferend, iname, varid
55  type(gd_nc_variable_entry):: ent
56  continue
57  stat = vtable_lookup(var, ent)
58  if (stat /= nf90_noerr) then
59  if (present(default)) value(:) = default
60  return
61  endif
62  ! 型と長さを取得
63  if (name(1:1) == gt_plus) then
64  iname = 2
65  varid = nf90_global
66  else
67  iname = 1
68  varid = ent%varid
69  endif
70  stat = nf90_inquire_attribute(ent%fileid, varid, &
71  name = name(iname:), xtype = xtype, len = attrlen)
72  if (stat /= nf90_noerr) then
73  if (present(default)) value(:) = default
74  return
75  endif
76  ! 文字型の場合は長さをコンマで分解した語数と読み替える
77  if (xtype == nf90_char) then
78  call get_attr(var, name, cbuffer, "", stat)
79  if (stat /= 0) return
80  call split(cbuffer, lbuffer, ", ")
81  attrlen = size(lbuffer)
82  endif
83  ! 結果を入れるところがなければ長さだけを伝え終了
84  if (size(value) == 0) then
85  if (xtype == nf90_char) deallocate(lbuffer)
86  stat = attrlen
87  return
88  endif
89  ! 型に応じて要求されただけ値を転送
90  xferend = min(size(value), attrlen)
91  if (present(default)) value(xferend+1: ) = default
92  if (xtype == nf90_char) then
93  do, i = 1, xferend
94  value(i) = stod(lbuffer(i))
95  enddo
96  deallocate(lbuffer)
97  stat = attrlen
98  return
99  else
100  allocate(rbuffer(attrlen), stat=stat)
101  if (stat /= 0) then
102  stat = nf90_enomem
103  return
104  endif
105  stat = nf90_get_att(ent%fileid, varid, name(iname:), rbuffer)
106  if (stat == nf90_noerr) then
107  value(1:xferend) = rbuffer(1:xferend)
108  stat = attrlen
109  endif
110  deallocate(rbuffer)
111  return
112  endif
integer function, public vtable_lookup(var, entry)
character, parameter, public gt_plus
Definition: dc_url.f90:92
文字型変数の操作.
Definition: dc_string.f90:24
種別型パラメタを提供します。
Definition: dc_types.f90:49
integer, parameter, public string
文字列を保持する 文字型変数の種別型パラメタ
Definition: dc_types.f90:118
Here is the call graph for this function:

◆ gdncattrgetreal()

subroutine gdncattrgetreal ( type(gd_nc_variable), intent(in)  var,
character(len = *), intent(in)  name,
real(sp), dimension(:), intent(out)  value,
integer, intent(out)  stat,
real(sp), intent(in), optional  default 
)

Definition at line 116 of file gdncattrgetnum.f90.

References dc_url::gt_plus, dc_types::sp, dc_types::string, and gtdata_netcdf_internal::vtable_lookup().

118  use netcdf, only: &
119  & nf90_noerr, &
120  & nf90_global, &
121  & nf90_char, &
122  & nf90_enomem, &
123  & nf90_inquire_attribute, &
124  & nf90_get_att
125  use dc_url, only: gt_plus
127  use dc_types, only: string, sp
128  use dc_string
129  implicit none
130  type(gd_nc_variable), intent(in) :: var
131  character(len = *), intent(in) :: name
132  integer, intent(out):: stat
133  real(SP), intent(out):: value(:)
134  real(SP), intent(in), optional:: default
135  real(SP), allocatable:: rbuffer(:)
136  character(STRING) :: cbuffer
137  character(STRING), pointer :: lbuffer(:)
138  integer :: attrlen, xtype, i, xferend, iname, varid
139  type(gd_nc_variable_entry):: ent
140  continue
141  stat = vtable_lookup(var, ent)
142  if (stat /= nf90_noerr) then
143  if (present(default)) value(:) = default
144  return
145  endif
146  ! 型と長さを取得
147  if (name(1:1) == gt_plus) then
148  iname = 2
149  varid = nf90_global
150  else
151  iname = 1
152  varid = ent%varid
153  endif
154  stat = nf90_inquire_attribute(ent%fileid, varid, &
155  name = name(iname:), xtype = xtype, len = attrlen)
156  if (stat /= nf90_noerr) then
157  if (present(default)) value(:) = default
158  return
159  endif
160  ! 文字型の場合は長さをコンマで分解した語数と読み替える
161  if (xtype == nf90_char) then
162  call get_attr(var, name, cbuffer, "", stat)
163  if (stat /= 0) return
164  call split(cbuffer, lbuffer, ", ")
165  attrlen = size(lbuffer)
166  endif
167  ! 結果を入れるところがなければ長さだけを伝え終了
168  if (size(value) == 0) then
169  if (xtype == nf90_char) deallocate(lbuffer)
170  stat = attrlen
171  return
172  endif
173  ! 型に応じて要求されただけ値を転送
174  xferend = min(size(value), attrlen)
175  if (present(default)) value(xferend+1: ) = default
176  if (xtype == nf90_char) then
177  do, i = 1, xferend
178  value(i) = stod(lbuffer(i))
179  enddo
180  deallocate(lbuffer)
181  stat = attrlen
182  return
183  else
184  allocate(rbuffer(attrlen), stat=stat)
185  if (stat /= 0) then
186  stat = nf90_enomem
187  return
188  endif
189  stat = nf90_get_att(ent%fileid, varid, name(iname:), rbuffer)
190  if (stat == nf90_noerr) then
191  value(1:xferend) = rbuffer(1:xferend)
192  stat = attrlen
193  endif
194  deallocate(rbuffer)
195  return
196  endif
integer function, public vtable_lookup(var, entry)
character, parameter, public gt_plus
Definition: dc_url.f90:92
文字型変数の操作.
Definition: dc_string.f90:24
種別型パラメタを提供します。
Definition: dc_types.f90:49
integer, parameter, public sp
単精度実数型変数
Definition: dc_types.f90:73
integer, parameter, public string
文字列を保持する 文字型変数の種別型パラメタ
Definition: dc_types.f90:118
Here is the call graph for this function: