historyautoputaxis.f90 File Reference

Go to the source code of this file.

Functions/Subroutines

subroutine historyautoputaxisreal (dim, array)
 
subroutine historyautoputaxisdouble (dim, array)
 
subroutine historyautoputaxisint (dim, array)
 

Function/Subroutine Documentation

◆ historyautoputaxisdouble()

subroutine historyautoputaxisdouble ( character(*), intent(in)  dim,
real(dp), dimension(:), intent(in)  array 
)

Definition at line 117 of file historyautoputaxis.f90.

References dc_trace::beginsub(), gtool_historyauto_internal::data_axes, dc_error::dc_enotinit, dc_error::dc_noerr, dc_types::dp, dc_trace::endsub(), dc_error::gt_eargsizemismatch, gtool_historyauto_internal::gthst_axes, dc_error::hst_enoaxisname, gtool_historyauto_internal::initialized, gtool_historyauto_internal::numdims, dc_error::storeerror(), dc_types::string, and dc_types::token.

117  !
118  ! 座標データを設定します.
119  !
120  ! Set data of an axis.
121  !
122 
126  use dc_trace, only: beginsub, endsub
127  use dc_error, only: storeerror, dc_noerr, &
129  use dc_types, only: dp, string, token
130 
131  implicit none
132  character(*), intent(in):: dim
133  ! 座標の名称.
134  !
135  ! ただし, ここで指定するもの
136  ! は, HistoryAutoCreate の *dims*
137  ! 既に指定されていなければなりません.
138  !
139  ! Name of axis.
140  !
141  ! Note that this value must be set
142  ! as "dims" of "HistoryAutoCreate".
143  !
144  real(DP), intent(in):: array(:)
145  ! 座標データ
146  !
147  ! データ型は整数, 単精度実数型,
148  ! 倍精度実数型のどれでもかまいません.
149  ! ただし, ファイルへ出力される際には,
150  ! HistoryAutoCreate の *xtypes* で指定した
151  ! データ型へ変換されます.
152  !
153  ! Data of axis
154  !
155  ! Integer, single or double precision are
156  ! acceptable as data type.
157  ! Note that when this is output to a file,
158  ! data type is converted into "xtypes"
159  ! specified in "HistoryAutoCreate"
160  !
161 
162  character(STRING):: name
163  integer:: stat, i
164  character(STRING):: cause_c
165  character(*), parameter:: subname = "HistoryAutoPutAxisDouble"
166  continue
167  call beginsub(subname, 'dim=<%c>', c1=trim(dim) )
168  stat = dc_noerr
169  cause_c = ""
170 
171  ! 初期設定チェック
172  ! Check initialization
173  !
174  if ( .not. initialized ) then
175  stat = dc_enotinit
176  cause_c = 'gtool_historyauto'
177  goto 999
178  end if
179 
180  do i = 1, numdims
181  call historyaxisinquire( &
182  & axis = gthst_axes(i), & ! (in)
183  & name = name ) ! (out)
184  if ( trim(dim) == trim(name) ) then
185  data_axes(i) % a_axis = array
186  goto 999
187  end if
188  end do
189 
190  stat = hst_enoaxisname
191  cause_c = dim
192 
193 999 continue
194  call storeerror(stat, subname, cause_c = cause_c)
195  call endsub(subname)
integer, parameter, public dc_enotinit
Definition: dc_error.f90:557
integer, parameter, public token
単語やキーワードを保持する文字型変数の種別型パラメタ
Definition: dc_types.f90:109
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
type(gt_history_axis_data), dimension(1:nf90_max_dims), target, save, public data_axes
type(gt_history_axis), dimension(1:nf90_max_dims), target, save, public gthst_axes
integer, parameter, public dp
倍精度実数型変数
Definition: dc_types.f90:83
subroutine, public beginsub(name, fmt, i, r, d, L, n, c1, c2, c3, ca, version)
Definition: dc_trace.f90:351
integer, parameter, public gt_eargsizemismatch
Definition: dc_error.f90:536
種別型パラメタを提供します。
Definition: dc_types.f90:49
integer, parameter, public hst_enoaxisname
Definition: dc_error.f90:589
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:

◆ historyautoputaxisint()

subroutine historyautoputaxisint ( character(*), intent(in)  dim,
integer, dimension(:), intent(in)  array 
)

Definition at line 203 of file historyautoputaxis.f90.

References dc_trace::beginsub(), gtool_historyauto_internal::data_axes, dc_error::dc_enotinit, dc_error::dc_noerr, dc_types::dp, dc_trace::endsub(), dc_error::gt_eargsizemismatch, gtool_historyauto_internal::gthst_axes, dc_error::hst_enoaxisname, gtool_historyauto_internal::initialized, gtool_historyauto_internal::numdims, dc_error::storeerror(), dc_types::string, and dc_types::token.

203  !
204  ! 座標データを設定します.
205  !
206  ! Set data of an axis.
207  !
208 
212  use dc_trace, only: beginsub, endsub
213  use dc_error, only: storeerror, dc_noerr, &
215  use dc_types, only: dp, string, token
216 
217  implicit none
218  character(*), intent(in):: dim
219  ! 座標の名称.
220  !
221  ! ただし, ここで指定するもの
222  ! は, HistoryAutoCreate の *dims*
223  ! 既に指定されていなければなりません.
224  !
225  ! Name of axis.
226  !
227  ! Note that this value must be set
228  ! as "dims" of "HistoryAutoCreate".
229  !
230  integer, intent(in):: array(:)
231  ! 座標データ
232  !
233  ! データ型は整数, 単精度実数型,
234  ! 倍精度実数型のどれでもかまいません.
235  ! ただし, ファイルへ出力される際には,
236  ! HistoryAutoCreate の *xtypes* で指定した
237  ! データ型へ変換されます.
238  !
239  ! Data of axis
240  !
241  ! Integer, single or double precision are
242  ! acceptable as data type.
243  ! Note that when this is output to a file,
244  ! data type is converted into "xtypes"
245  ! specified in "HistoryAutoCreate"
246  !
247 
248  character(STRING):: name
249  integer:: stat, i
250  character(STRING):: cause_c
251  character(*), parameter:: subname = "HistoryAutoPutAxisInt"
252  continue
253  call beginsub(subname, 'dim=<%c>', c1=trim(dim) )
254  stat = dc_noerr
255  cause_c = ""
256 
257  ! 初期設定チェック
258  ! Check initialization
259  !
260  if ( .not. initialized ) then
261  stat = dc_enotinit
262  cause_c = 'gtool_historyauto'
263  goto 999
264  end if
265 
266  do i = 1, numdims
267  call historyaxisinquire( &
268  & axis = gthst_axes(i), & ! (in)
269  & name = name ) ! (out)
270  if ( trim(dim) == trim(name) ) then
271  data_axes(i) % a_axis = array
272  goto 999
273  end if
274  end do
275 
276  stat = hst_enoaxisname
277  cause_c = dim
278 
279 999 continue
280  call storeerror(stat, subname, cause_c = cause_c)
281  call endsub(subname)
integer, parameter, public dc_enotinit
Definition: dc_error.f90:557
integer, parameter, public token
単語やキーワードを保持する文字型変数の種別型パラメタ
Definition: dc_types.f90:109
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
type(gt_history_axis_data), dimension(1:nf90_max_dims), target, save, public data_axes
type(gt_history_axis), dimension(1:nf90_max_dims), target, save, public gthst_axes
integer, parameter, public dp
倍精度実数型変数
Definition: dc_types.f90:83
subroutine, public beginsub(name, fmt, i, r, d, L, n, c1, c2, c3, ca, version)
Definition: dc_trace.f90:351
integer, parameter, public gt_eargsizemismatch
Definition: dc_error.f90:536
種別型パラメタを提供します。
Definition: dc_types.f90:49
integer, parameter, public hst_enoaxisname
Definition: dc_error.f90:589
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:

◆ historyautoputaxisreal()

subroutine historyautoputaxisreal ( character(*), intent(in)  dim,
real, dimension(:), intent(in)  array 
)

Definition at line 31 of file historyautoputaxis.f90.

References dc_trace::beginsub(), gtool_historyauto_internal::data_axes, dc_error::dc_enotinit, dc_error::dc_noerr, dc_types::dp, dc_trace::endsub(), dc_error::gt_eargsizemismatch, gtool_historyauto_internal::gthst_axes, dc_error::hst_enoaxisname, gtool_historyauto_internal::initialized, gtool_historyauto_internal::numdims, dc_error::storeerror(), dc_types::string, and dc_types::token.

31  !
32  ! 座標データを設定します.
33  !
34  ! Set data of an axis.
35  !
36 
40  use dc_trace, only: beginsub, endsub
41  use dc_error, only: storeerror, dc_noerr, &
43  use dc_types, only: dp, string, token
44 
45  implicit none
46  character(*), intent(in):: dim
47  ! 座標の名称.
48  !
49  ! ただし, ここで指定するもの
50  ! は, HistoryAutoCreate の *dims*
51  ! 既に指定されていなければなりません.
52  !
53  ! Name of axis.
54  !
55  ! Note that this value must be set
56  ! as "dims" of "HistoryAutoCreate".
57  !
58  real, intent(in):: array(:)
59  ! 座標データ
60  !
61  ! データ型は整数, 単精度実数型,
62  ! 倍精度実数型のどれでもかまいません.
63  ! ただし, ファイルへ出力される際には,
64  ! HistoryAutoCreate の *xtypes* で指定した
65  ! データ型へ変換されます.
66  !
67  ! Data of axis
68  !
69  ! Integer, single or double precision are
70  ! acceptable as data type.
71  ! Note that when this is output to a file,
72  ! data type is converted into "xtypes"
73  ! specified in "HistoryAutoCreate"
74  !
75 
76  character(STRING):: name
77  integer:: stat, i
78  character(STRING):: cause_c
79  character(*), parameter:: subname = "HistoryAutoPutAxisReal"
80  continue
81  call beginsub(subname, 'dim=<%c>', c1=trim(dim) )
82  stat = dc_noerr
83  cause_c = ""
84 
85  ! 初期設定チェック
86  ! Check initialization
87  !
88  if ( .not. initialized ) then
89  stat = dc_enotinit
90  cause_c = 'gtool_historyauto'
91  goto 999
92  end if
93 
94  do i = 1, numdims
95  call historyaxisinquire( &
96  & axis = gthst_axes(i), & ! (in)
97  & name = name ) ! (out)
98  if ( trim(dim) == trim(name) ) then
99  data_axes(i) % a_axis = array
100  goto 999
101  end if
102  end do
103 
104  stat = hst_enoaxisname
105  cause_c = dim
106 
107 999 continue
108  call storeerror(stat, subname, cause_c = cause_c)
109  call endsub(subname)
integer, parameter, public dc_enotinit
Definition: dc_error.f90:557
integer, parameter, public token
単語やキーワードを保持する文字型変数の種別型パラメタ
Definition: dc_types.f90:109
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
type(gt_history_axis_data), dimension(1:nf90_max_dims), target, save, public data_axes
type(gt_history_axis), dimension(1:nf90_max_dims), target, save, public gthst_axes
integer, parameter, public dp
倍精度実数型変数
Definition: dc_types.f90:83
subroutine, public beginsub(name, fmt, i, r, d, L, n, c1, c2, c3, ca, version)
Definition: dc_trace.f90:351
integer, parameter, public gt_eargsizemismatch
Definition: dc_error.f90:536
種別型パラメタを提供します。
Definition: dc_types.f90:49
integer, parameter, public hst_enoaxisname
Definition: dc_error.f90:589
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: