gtvargetnum.f90
Go to the documentation of this file.
1 ! -*- coding: utf-8; mode: f90 -*-
2 !-------------------------------------------------------------------------------------
3 ! Copyright (c) 2000-2016 Gtool Development Group. All rights reserved.
4 !-------------------------------------------------------------------------------------
5 ! ** Important**
6 !
7 ! This file is generated from gtvargetnum.erb by ERB included Ruby 2.3.1.
8 ! Please do not edit this file directly. @see "gtvargetnum.erb"
9 !-------------------------------------------------------------------------------------
10 !
11 != 固定長配列への数値データの入力
12 !
13 ! Authors:: Yasuhiro MORIKAWA, Eizi TOYODA
14 ! Version:: $Id: gtvargetnum.rb2f90,v 1.5 2009-05-25 09:55:58 morikawa Exp $
15 ! Tag Name:: $Name: $
16 ! Copyright:: Copyright (C) GFD Dennou Club, 2000-2005. All rights reserved.
17 ! License:: See COPYRIGHT[link:../../COPYRIGHT]
18 !
19 ! 以下のサブルーチン、関数は gtdata_generic から gtdata_generic#Get
20 ! として提供されます。
21 !
22 
23 !== 固定長配列への数値データの入力
24 !
25 ! 変数 *var* から *value* に数値データが入力されます。
26 ! *nvalue* には配列長を代入する必要があります。
27 !
28 ! 数値データ入力の際にエラーが生じた場合、メッセージを出力
29 ! してプログラムは強制終了します。*err* を与えてある場合には
30 ! の引数に .true. が返り、プログラムは終了しません。
31 !
32 ! 入力しようとするデータの型が引数の型と異なる場合、データは引数の
33 ! 型に変換されます。 この変換は netCDF の機能を用いています。
34 ! 詳しくは {netCDF 日本語版マニュアル}[link:../xref.htm#label-10]
35 ! の 3.3 型変換 を参照してください。
36 !
37 ! *Get* は複数のサブルーチンの総称名であり、
38 ! *value* にポインタ型の配列を与えることも可能です。上記の
39 ! サブルーチンを参照してください。
40 
41 subroutine gtvargetdouble(var, value, nvalue, err)
43  use gtdata_internal_map, only: &
44  & var_class, &
45  & vtb_class_netcdf, &
47  & gtvar_dump
48  use gtdata_netcdf_generic, only: get
50  use dc_types, only: string, dp
52  implicit none
53  type(gt_variable), intent(in) :: var
54  real(DP), intent(out) :: value(*)
55  integer, intent(in) :: nvalue
56  logical, intent(out), optional :: err
57  integer :: class, cid, stat
58  integer , pointer :: specs(:, :)
59  character(STRING):: cause_c
60  character(len = *), parameter:: subname = 'GTVarGetDouble'
61  continue
62  call var_class(var, class, cid)
63  stat = dc_noerr
64  cause_c = ''
65  if (nvalue < 1) then
66  stat = dc_enegative
67  cause_c = 'nvalue'
68  goto 999
69  end if
70  call gtvar_dump(var)
71  call map_to_internal_specs(var, specs)
72  if (class == vtb_class_netcdf) then
73  call get(gd_nc_variable(cid), start=specs(:, 1), count=specs(:, 2), &
74  & stride=specs(:, 3), imap=specs(:, 4), siz=nvalue, value=value, &
75  & iostat=stat)
76  else
77  stat = gt_efake
78  endif
79  if (associated(specs)) deallocate(specs)
80 999 continue
81  call storeerror(stat, subname, cause_c = cause_c, err = err)
82 end subroutine gtvargetdouble
83 
84 subroutine gtvargetreal(var, value, nvalue, err)
86  use gtdata_internal_map, only: &
87  & var_class, &
88  & vtb_class_netcdf, &
90  & gtvar_dump
91  use gtdata_netcdf_generic, only: get
93  use dc_types, only: string ,sp
95  implicit none
96  type(gt_variable), intent(in) :: var
97  real(SP), intent(out) :: value(*)
98  integer, intent(in) :: nvalue
99  logical, intent(out), optional :: err
100  integer :: class, cid, stat
101  integer , pointer :: specs(:, :)
102  character(STRING):: cause_c
103  character(len = *), parameter:: subname = 'GTVarGetDouble'
104  continue
105  call var_class(var, class, cid)
106  stat = dc_noerr
107  cause_c = ''
108  if (nvalue < 1) then
109  stat = dc_enegative
110  cause_c = 'nvalue'
111  goto 999
112  end if
113  call gtvar_dump(var)
114  call map_to_internal_specs(var, specs)
115  if (class == vtb_class_netcdf) then
116  call get(gd_nc_variable(cid), start=specs(:, 1), count=specs(:, 2), &
117  & stride=specs(:, 3), imap=specs(:, 4), siz=nvalue, value=value, &
118  & iostat=stat)
119  else
120  stat = gt_efake
121  endif
122  if (associated(specs)) deallocate(specs)
123 999 continue
124  call storeerror(stat, subname, cause_c = cause_c, err = err)
125 end subroutine gtvargetreal
126 
127 subroutine gtvargetint(var, value, nvalue, err)
129  use gtdata_internal_map, only: &
130  & var_class, &
131  & vtb_class_netcdf, &
133  & gtvar_dump
134  use gtdata_netcdf_generic, only: get
136  use dc_types, only: string
138  implicit none
139  type(gt_variable), intent(in) :: var
140  integer, intent(out) :: value(*)
141  integer, intent(in) :: nvalue
142  logical, intent(out), optional :: err
143  integer :: class, cid, stat
144  integer , pointer :: specs(:, :)
145  character(STRING):: cause_c
146  character(len = *), parameter:: subname = 'GTVarGetDouble'
147  continue
148  call var_class(var, class, cid)
149  stat = dc_noerr
150  cause_c = ''
151  if (nvalue < 1) then
152  stat = dc_enegative
153  cause_c = 'nvalue'
154  goto 999
155  end if
156  call gtvar_dump(var)
157  call map_to_internal_specs(var, specs)
158  if (class == vtb_class_netcdf) then
159  call get(gd_nc_variable(cid), start=specs(:, 1), count=specs(:, 2), &
160  & stride=specs(:, 3), imap=specs(:, 4), siz=nvalue, value=value, &
161  & iostat=stat)
162  else
163  stat = gt_efake
164  endif
165  if (associated(specs)) deallocate(specs)
166 999 continue
167  call storeerror(stat, subname, cause_c = cause_c, err = err)
168 end subroutine gtvargetint
169 
subroutine, public map_to_internal_specs(var, specs, ndims)
integer, parameter, public gt_efake
Definition: dc_error.f90:523
subroutine, public storeerror(number, where, err, cause_c, cause_i)
Definition: dc_error.f90:830
integer, parameter, public dc_noerr
Definition: dc_error.f90:509
subroutine gtvargetdouble(var, value, nvalue, err)
Definition: gtvargetnum.f90:42
integer, parameter, public dp
倍精度実数型変数
Definition: dc_types.f90:83
subroutine gtvargetint(var, value, nvalue, err)
種別型パラメタを提供します。
Definition: dc_types.f90:49
integer, parameter, public sp
単精度実数型変数
Definition: dc_types.f90:73
subroutine, public var_class(var, class, cid)
subroutine gtvargetreal(var, value, nvalue, err)
Definition: gtvargetnum.f90:85
integer, parameter, public dc_enegative
Definition: dc_error.f90:568
integer, parameter, public string
文字列を保持する 文字型変数の種別型パラメタ
Definition: dc_types.f90:118