historyautoaddweight.f90 File Reference

Go to the source code of this file.

Functions/Subroutines

subroutine historyautoaddweightreal (dim, weight, units, xtype)
 
subroutine historyautoaddweightdouble (dim, weight, units, xtype)
 
subroutine historyautoaddweightint (dim, weight, units, xtype)
 

Function/Subroutine Documentation

◆ historyautoaddweightdouble()

subroutine historyautoaddweightdouble ( character(*), intent(in)  dim,
real(dp), dimension(:), intent(in)  weight,
character(*), intent(in), optional  units,
character(*), intent(in), optional  xtype 
)

Definition at line 164 of file historyautoaddweight.f90.

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

164  !
165  ! 座標の重みデータを設定します.
166  !
167  ! Set weights of axes.
168  !
169 
174  use dc_trace, only: beginsub, endsub
175  use dc_error, only: storeerror, dc_noerr, &
177  use dc_types, only: dp, string, token
178 
179  implicit none
180  character(*), intent(in):: dim
181  ! 座標重みを設定する座標の名称.
182  !
183  ! ただし, ここで指定するもの
184  ! は, HistoryAutoCreate の *dims*
185  ! 既に指定されていなければなりません.
186  !
187  ! Name of axis to which "weight" are set.
188  !
189  ! Note that this value must be set
190  ! as "dims" of "HistoryAutoCreate".
191  !
192  real(DP), intent(in):: weight(:)
193  ! 座標重みデータ.
194  !
195  ! データ型は整数, 単精度実数型,
196  ! 倍精度実数型のどれでもかまいません.
197  ! ただし, ファイルへ出力される際には,
198  ! xtype もしくは座標データの型へと
199  ! 変換されます.
200  !
201  ! Weight of axis.
202  !
203  ! Integer, single or double precision are
204  ! acceptable as data type.
205  ! Note that when this is output to a file,
206  ! data type is converted into "xtype" or
207  ! type of the axis.
208  !
209  character(*), intent(in), optional:: units
210  ! 座標重みの単位.
211  ! 省略した場合には, 座標の単位が
212  ! 使用されます.
213  !
214  ! Units of axis weight.
215  ! If this argument is omitted,
216  ! unit of the dimension is used.
217  !
218  character(*), intent(in), optional:: xtype
219  ! 座標重みのデータ型.
220  ! 省略した場合には, 座標のデータ型が
221  ! 使用されます.
222  !
223  ! Data type of weight of the dimension.
224  ! If this argument is omitted,
225  ! data type of the dimension is used.
226  !
227 
228  character(STRING):: name, longname
229  character(TOKEN):: dim_units, dim_xtype
230  integer:: dim_size
231  integer:: stat, i
232  character(STRING):: cause_c
233  character(*), parameter:: subname = "HistoryAutoAddWeightDouble"
234  continue
235  call beginsub(subname, 'dim=<%c>', c1=trim(dim) )
236  stat = dc_noerr
237  cause_c = ""
238 
239  ! 初期設定チェック
240  ! Check initialization
241  !
242  if ( .not. initialized ) then
243  stat = dc_enotinit
244  cause_c = 'gtool_historyauto'
245  goto 999
246  end if
247 
248  do i = 1, numdims
249  call historyaxisinquire( &
250  & axis = gthst_axes(i), & ! (in)
251  & name = name, & ! (out)
252  & size = dim_size, & ! (out)
253  & longname = longname, & ! (out)
254  & units = dim_units, & ! (out)
255  & xtype = dim_xtype ) ! (out)
256  if ( trim(dim) == trim(name) ) then
257  if ( dim_size /= size(weight) ) then
258  stat = gt_eargsizemismatch
259  cause_c = 'weight'
260  end if
261  if ( present(units) ) dim_units = units
262  if ( present(xtype) ) dim_xtype = xtype
263  call historyvarinfocreate( &
264  & varinfo = gthst_weights(numwgts + 1), & ! (out)
265  & name = trim(dim) // wgtsuf, dims = (/ dim /), & ! (in)
266  & longname = 'weight for integration or average in ' // &
267  & trim(longname), & ! (in)
268  & units = dim_units, xtype = dim_xtype ) ! (in)
269 
270  call historyaxisaddattr( &
271  & axis = gthst_axes(i), & ! (inout)
272  & attrname = 'gt_calc_weight', & ! (in)
273  & value = trim(dim) // wgtsuf ) ! (in)
274 
275  allocate( data_weights(numwgts + 1) % a_axis( dim_size ) )
276  data_weights(numwgts + 1) % a_axis = weight
277 
278  numwgts = numwgts + 1
279  goto 999
280  end if
281  end do
282 
283  stat = hst_enoaxisname
284  cause_c = dim
285 
286 999 continue
287  call storeerror(stat, subname, cause_c = cause_c)
288  call endsub(subname)
integer, parameter, public dc_enotinit
Definition: dc_error.f90:557
character(*), parameter, public wgtsuf
integer, parameter, public token
単語やキーワードを保持する文字型変数の種別型パラメタ
Definition: dc_types.f90:109
type(gt_history_axis_data), dimension(1:nf90_max_dims), target, save, public data_weights
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_varinfo), dimension(1:nf90_max_dims), save, public gthst_weights
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:

◆ historyautoaddweightint()

subroutine historyautoaddweightint ( character(*), intent(in)  dim,
integer, dimension(:), intent(in)  weight,
character(*), intent(in), optional  units,
character(*), intent(in), optional  xtype 
)

Definition at line 296 of file historyautoaddweight.f90.

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

296  !
297  ! 座標の重みデータを設定します.
298  !
299  ! Set weights of axes.
300  !
301 
306  use dc_trace, only: beginsub, endsub
307  use dc_error, only: storeerror, dc_noerr, &
309  use dc_types, only: dp, string, token
310 
311  implicit none
312  character(*), intent(in):: dim
313  ! 座標重みを設定する座標の名称.
314  !
315  ! ただし, ここで指定するもの
316  ! は, HistoryAutoCreate の *dims*
317  ! 既に指定されていなければなりません.
318  !
319  ! Name of axis to which "weight" are set.
320  !
321  ! Note that this value must be set
322  ! as "dims" of "HistoryAutoCreate".
323  !
324  integer, intent(in):: weight(:)
325  ! 座標重みデータ.
326  !
327  ! データ型は整数, 単精度実数型,
328  ! 倍精度実数型のどれでもかまいません.
329  ! ただし, ファイルへ出力される際には,
330  ! xtype もしくは座標データの型へと
331  ! 変換されます.
332  !
333  ! Weight of axis.
334  !
335  ! Integer, single or double precision are
336  ! acceptable as data type.
337  ! Note that when this is output to a file,
338  ! data type is converted into "xtype" or
339  ! type of the axis.
340  !
341  character(*), intent(in), optional:: units
342  ! 座標重みの単位.
343  ! 省略した場合には, 座標の単位が
344  ! 使用されます.
345  !
346  ! Units of axis weight.
347  ! If this argument is omitted,
348  ! unit of the dimension is used.
349  !
350  character(*), intent(in), optional:: xtype
351  ! 座標重みのデータ型.
352  ! 省略した場合には, 座標のデータ型が
353  ! 使用されます.
354  !
355  ! Data type of weight of the dimension.
356  ! If this argument is omitted,
357  ! data type of the dimension is used.
358  !
359 
360  character(STRING):: name, longname
361  character(TOKEN):: dim_units, dim_xtype
362  integer:: dim_size
363  integer:: stat, i
364  character(STRING):: cause_c
365  character(*), parameter:: subname = "HistoryAutoAddWeightInt"
366  continue
367  call beginsub(subname, 'dim=<%c>', c1=trim(dim) )
368  stat = dc_noerr
369  cause_c = ""
370 
371  ! 初期設定チェック
372  ! Check initialization
373  !
374  if ( .not. initialized ) then
375  stat = dc_enotinit
376  cause_c = 'gtool_historyauto'
377  goto 999
378  end if
379 
380  do i = 1, numdims
381  call historyaxisinquire( &
382  & axis = gthst_axes(i), & ! (in)
383  & name = name, & ! (out)
384  & size = dim_size, & ! (out)
385  & longname = longname, & ! (out)
386  & units = dim_units, & ! (out)
387  & xtype = dim_xtype ) ! (out)
388  if ( trim(dim) == trim(name) ) then
389  if ( dim_size /= size(weight) ) then
390  stat = gt_eargsizemismatch
391  cause_c = 'weight'
392  end if
393  if ( present(units) ) dim_units = units
394  if ( present(xtype) ) dim_xtype = xtype
395  call historyvarinfocreate( &
396  & varinfo = gthst_weights(numwgts + 1), & ! (out)
397  & name = trim(dim) // wgtsuf, dims = (/ dim /), & ! (in)
398  & longname = 'weight for integration or average in ' // &
399  & trim(longname), & ! (in)
400  & units = dim_units, xtype = dim_xtype ) ! (in)
401 
402  call historyaxisaddattr( &
403  & axis = gthst_axes(i), & ! (inout)
404  & attrname = 'gt_calc_weight', & ! (in)
405  & value = trim(dim) // wgtsuf ) ! (in)
406 
407  allocate( data_weights(numwgts + 1) % a_axis( dim_size ) )
408  data_weights(numwgts + 1) % a_axis = weight
409 
410  numwgts = numwgts + 1
411  goto 999
412  end if
413  end do
414 
415  stat = hst_enoaxisname
416  cause_c = dim
417 
418 999 continue
419  call storeerror(stat, subname, cause_c = cause_c)
420  call endsub(subname)
integer, parameter, public dc_enotinit
Definition: dc_error.f90:557
character(*), parameter, public wgtsuf
integer, parameter, public token
単語やキーワードを保持する文字型変数の種別型パラメタ
Definition: dc_types.f90:109
type(gt_history_axis_data), dimension(1:nf90_max_dims), target, save, public data_weights
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_varinfo), dimension(1:nf90_max_dims), save, public gthst_weights
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:

◆ historyautoaddweightreal()

subroutine historyautoaddweightreal ( character(*), intent(in)  dim,
real, dimension(:), intent(in)  weight,
character(*), intent(in), optional  units,
character(*), intent(in), optional  xtype 
)

Definition at line 32 of file historyautoaddweight.f90.

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

32  !
33  ! 座標の重みデータを設定します.
34  !
35  ! Set weights of axes.
36  !
37 
42  use dc_trace, only: beginsub, endsub
43  use dc_error, only: storeerror, dc_noerr, &
45  use dc_types, only: dp, string, token
46 
47  implicit none
48  character(*), intent(in):: dim
49  ! 座標重みを設定する座標の名称.
50  !
51  ! ただし, ここで指定するもの
52  ! は, HistoryAutoCreate の *dims*
53  ! 既に指定されていなければなりません.
54  !
55  ! Name of axis to which "weight" are set.
56  !
57  ! Note that this value must be set
58  ! as "dims" of "HistoryAutoCreate".
59  !
60  real, intent(in):: weight(:)
61  ! 座標重みデータ.
62  !
63  ! データ型は整数, 単精度実数型,
64  ! 倍精度実数型のどれでもかまいません.
65  ! ただし, ファイルへ出力される際には,
66  ! xtype もしくは座標データの型へと
67  ! 変換されます.
68  !
69  ! Weight of axis.
70  !
71  ! Integer, single or double precision are
72  ! acceptable as data type.
73  ! Note that when this is output to a file,
74  ! data type is converted into "xtype" or
75  ! type of the axis.
76  !
77  character(*), intent(in), optional:: units
78  ! 座標重みの単位.
79  ! 省略した場合には, 座標の単位が
80  ! 使用されます.
81  !
82  ! Units of axis weight.
83  ! If this argument is omitted,
84  ! unit of the dimension is used.
85  !
86  character(*), intent(in), optional:: xtype
87  ! 座標重みのデータ型.
88  ! 省略した場合には, 座標のデータ型が
89  ! 使用されます.
90  !
91  ! Data type of weight of the dimension.
92  ! If this argument is omitted,
93  ! data type of the dimension is used.
94  !
95 
96  character(STRING):: name, longname
97  character(TOKEN):: dim_units, dim_xtype
98  integer:: dim_size
99  integer:: stat, i
100  character(STRING):: cause_c
101  character(*), parameter:: subname = "HistoryAutoAddWeightReal"
102  continue
103  call beginsub(subname, 'dim=<%c>', c1=trim(dim) )
104  stat = dc_noerr
105  cause_c = ""
106 
107  ! 初期設定チェック
108  ! Check initialization
109  !
110  if ( .not. initialized ) then
111  stat = dc_enotinit
112  cause_c = 'gtool_historyauto'
113  goto 999
114  end if
115 
116  do i = 1, numdims
117  call historyaxisinquire( &
118  & axis = gthst_axes(i), & ! (in)
119  & name = name, & ! (out)
120  & size = dim_size, & ! (out)
121  & longname = longname, & ! (out)
122  & units = dim_units, & ! (out)
123  & xtype = dim_xtype ) ! (out)
124  if ( trim(dim) == trim(name) ) then
125  if ( dim_size /= size(weight) ) then
126  stat = gt_eargsizemismatch
127  cause_c = 'weight'
128  end if
129  if ( present(units) ) dim_units = units
130  if ( present(xtype) ) dim_xtype = xtype
131  call historyvarinfocreate( &
132  & varinfo = gthst_weights(numwgts + 1), & ! (out)
133  & name = trim(dim) // wgtsuf, dims = (/ dim /), & ! (in)
134  & longname = 'weight for integration or average in ' // &
135  & trim(longname), & ! (in)
136  & units = dim_units, xtype = dim_xtype ) ! (in)
137 
138  call historyaxisaddattr( &
139  & axis = gthst_axes(i), & ! (inout)
140  & attrname = 'gt_calc_weight', & ! (in)
141  & value = trim(dim) // wgtsuf ) ! (in)
142 
143  allocate( data_weights(numwgts + 1) % a_axis( dim_size ) )
144  data_weights(numwgts + 1) % a_axis = weight
145 
146  numwgts = numwgts + 1
147  goto 999
148  end if
149  end do
150 
151  stat = hst_enoaxisname
152  cause_c = dim
153 
154 999 continue
155  call storeerror(stat, subname, cause_c = cause_c)
156  call endsub(subname)
integer, parameter, public dc_enotinit
Definition: dc_error.f90:557
character(*), parameter, public wgtsuf
integer, parameter, public token
単語やキーワードを保持する文字型変数の種別型パラメタ
Definition: dc_types.f90:109
type(gt_history_axis_data), dimension(1:nf90_max_dims), target, save, public data_weights
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_varinfo), dimension(1:nf90_max_dims), save, public gthst_weights
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: