gtvarlimit.f90 File Reference

Go to the source code of this file.

Functions/Subroutines

subroutine gtvarlimit_iiii (var, dimord, start, count, stride, err)
 
subroutine gtvarlimit (var, string, err)
 
subroutine limit_one (string)
 
subroutine region_spec (dimord, string, start, count, stride)
 
real function value_to_index (dimord, value)
 

Function/Subroutine Documentation

◆ gtvarlimit()

subroutine gtvarlimit ( type(gt_variable), intent(inout)  var,
character(len = *), intent(in)  string,
logical, intent(out), optional  err 
)

Definition at line 147 of file gtvarlimit.f90.

References dc_trace::beginsub(), dc_trace::endsub(), dc_url::gt_comma, gtdata_internal_map::gtvar_dump(), and limit_one().

147  !
148  !== 入出力範囲の拘束 (文字列で指定)
149  !
150  ! 変数 *var* 次元の入出力範囲を拘束します。
151  ! *Limit* は 2 つのサブルーチンの総称名であり、
152  ! 別の指定方法もあります。まずは上記のサブルーチンを参照してください。
153  !
154  ! 指定方法は、*string* に
155  ! {gtool4 netCDF 規約}[link:../xref.htm#label-6] の
156  ! 「5.4 コンマ記法」に述べられる範囲指定表現を用います。
157  ! 凡例を以下に挙げます。
158  !
159  ! <dim>=<lower>
160  !
161  ! <dim>=<lower>:<upper>
162  !
163  ! <dim>=<lower>:<upper>:<stride>
164  !
165  ! ここで、<tt><dim></tt> は次元番号または次元名であり、
166  ! <tt><lower></tt>, <tt><upper></tt>
167  ! は座標値または "<tt>^</tt>" を前置した格子番号です。
168  ! <tt><stride></tt> は格子数です。
169  !
170  ! エラーが生じた場合、メッセージを出力
171  ! してプログラムは強制終了します。*err* を与えてある場合には
172  ! の引数に .true. が返り、プログラムは終了しません。
173  !
174  use gtdata_types, only: gt_variable
175  use dc_trace, only: beginsub, endsub
176  use dc_url, only: gt_comma
178  type(gt_variable), intent(inout):: var
179  character(len = *), intent(in) :: string
180  logical, intent(out), optional :: err
181  integer:: is, ie
182 continue
183  call beginsub('GTVarLimit', 'var=%d lim=<%c>', i=(/var%mapid/), c1=trim(string))
184  call gtvar_dump(var)
185  ! コンマで区切って解釈
186  is = 1
187  do
188  ie = index(string(is: ), gt_comma)
189  if (ie == 0) exit
190  call limit_one(string(is: is+ie-2))
191  is = is + ie
192  if (is > len(string)) exit
193  enddo
194  call limit_one(string(is: ))
195  if (present(err)) err = .false.
196  call endsub('GTVarLimit')
197  return
198 contains
199 
200  subroutine limit_one(string)
201  use dc_url, only: gt_equal
202  use dc_string, only: strieq, stoi
205  character(len = *), intent(in):: string
206  integer:: equal, dimord
207  integer:: start, count, stride, strhead
208  logical:: myerr
209 
210  if (string == '') return
211 
212  strhead = 4
213  if (len(string) < 4) strhead = len(string)
214 
215  if (strieq(string(1:strhead), "IGN:")) then
216  ! 隠蔽型指定子 ign:<dim> または ign:<dim>=<start>
217  equal = index(string, gt_equal)
218  if (equal == 0) then
219  start = 1
220  else
221  start = stoi(string(equal+1: ), default=1)
222  endif
223  dimord = dimname_to_dimord(var, string(5: equal-1))
224  call limit(var, dimord, start, 1, 1, err)
225  call del_dim(var, dimord, myerr)
226  return
227  endif
228 
229  ! 限定型指定子 <dim>=<start>:<finish>:<stride>
230  ! いまは実装がバグっていて <start>:<count>:<stride> になってる
231  !
232  equal = index(string, gt_equal)
233  if (equal == 0) return
234  dimord = dimname_to_dimord(var, string(1: equal-1))
235  if (dimord <= 0) return
236  !
237  call region_spec(dimord, string(equal+1: ), start, count, stride)
238  call limit(var, dimord, start, count, stride, err)
239  end subroutine limit_one
240 
241  !
242  ! 範囲指定の = のあとを : で区切ってマップにいれる
243  !
244  subroutine region_spec(dimord, string, start, count, stride)
245  use dc_types, only: token
246  use dc_string, only: index_ofs, stoi
247  use dc_url, only: gt_circumflex, gt_colon
248  use gtdata_internal_map, only: dimrange
249  integer, intent(in):: dimord
250  integer, intent(out):: start, count, stride
251  character(len = *), intent(in):: string
252  integer:: colon, prev_colon, finish, dimlo, dimhi
253  character(len = token):: val(3)
254  continue
255  colon = index(string, gt_colon)
256  if (colon == 0) then
257  ! コロンがない場合は上下端に同じ値
258  val(1) = string(1: )
259  val(2) = val(1)
260  val(3) = ""
261  else
262  val(1) = string(1: colon - 1)
263  prev_colon = colon
264  colon = index_ofs(string, colon + 1, gt_colon)
265  if (colon > 0) then
266  val(2) = string(prev_colon + 1: colon - 1)
267  val(3) = string(colon + 1: )
268  else
269  val(2) = string(prev_colon + 1: )
270  val(3) = ""
271  endif
272  endif
273  if (val(3) == "") val(3) = "^1"
274 
275  if (val(1)(1:1) == gt_circumflex) then
276  start = stoi(val(1)(2: ))
277  else if (val(1) == val(2)) then
278  start = nint(value_to_index(dimord, val(1)))
279  else
280  start = floor(value_to_index(dimord, val(1)))
281  endif
282  if (val(2) == val(1)) then
283  finish = start
284  else if (val(2)(1:1) == gt_circumflex) then
285  finish = stoi(val(2)(2: ))
286  else
287  finish = ceiling(value_to_index(dimord, val(2)))
288  endif
289 
290  call dimrange(var, dimord, dimlo, dimhi)
291  start = min(max(dimlo, start), dimhi)
292  finish = min(max(dimlo, finish), dimhi)
293  count = abs(finish - start) + 1
294 
295  if (val(3)(1:1) == gt_circumflex) then
296  stride = stoi(val(3)(2: ))
297  else
298  stride = stoi(val(3))
299  endif
300  stride = sign(stride, finish - start)
301  end subroutine region_spec
302 
303  real function value_to_index(dimord, value) result(result)
304  !
305  ! GTVarLimit の引数 *var* に格納される変数の次元 *dimord*
306  ! に格納されるデータのうち, *value* が格納される
307  ! 格子番号を整数値にして返します.
308  !
309  ! 例えば次元に以下のデータが格納されているとします.
310  !
311  ! 0.05 0.1 0.15 0.20 0.25 0.30
312  !
313  ! この場合, *value* に 0.15 が与えられれば戻り値は 3. となります.
314  ! また *value* に 0.225 が与えられれば戻り値は 4.5 となります.
315  !
316  !
317  use gtdata_types, only: gt_variable
318  use gtdata_generic, only: get, open, close
319  use dc_string, only: stod
320  use dc_trace, only: beginsub, endsub, dbgmessage
321  integer, intent(in):: dimord
322  character(len = *), intent(in):: value
323  type(gt_variable):: axisvar
324  real, pointer:: axisval(:) => null()
325  real:: val
326  integer:: i
327  continue
328 
329  call beginsub('value_to_index', 'var=%d dimord=%d value=%c', &
330  & i=(/var%mapid, dimord/), c1=trim(value))
331 
332  call open(axisvar, var, dimord, count_compact=.true.)
333  nullify(axisval)
334  call get(axisvar, axisval)
335  call close(axisvar)
336  if (.not. associated(axisval)) then
337  result = -1.0
338  return
339  else if (size(axisval) < 2) then
340  result = 1.0
341  goto 900
342  endif
343 
344  val = stod(value)
345 
346  ! call DbgMessage('value=%f axis=(/%*r/)', r=(/val, axisval(:)/), &
347  ! & n=(/size(axisval)/))
348 
349  do, i = 1, size(axisval) - 1
350  if (axisval(i + 1) == axisval(i)) then
351  result = real(i) + 0.5
352  goto 900
353  endif
354  result = i + (val - axisval(i)) / (axisval(i + 1) - axisval(i))
355  if (result <= (i + 1)) goto 900
356  enddo
357 
358 900 continue
359  call endsub('value_to_index', '(%c) = %r', &
360  & c1=trim(value), r=(/result/))
361  deallocate(axisval)
362  end function value_to_index
363 
character, parameter, public gt_comma
Definition: dc_url.f90:85
integer, parameter, public token
単語やキーワードを保持する文字型変数の種別型パラメタ
Definition: dc_types.f90:109
character, parameter, public gt_equal
Definition: dc_url.f90:87
subroutine limit_one(string)
Definition: gtvarlimit.f90:201
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
character, parameter, public gt_circumflex
Definition: dc_url.f90:89
character, parameter, public gt_colon
Definition: dc_url.f90:83
文字型変数の操作.
Definition: dc_string.f90:24
種別型パラメタを提供します。
Definition: dc_types.f90:49
subroutine region_spec(dimord, string, start, count, stride)
Definition: gtvarlimit.f90:245
real function value_to_index(dimord, value)
Definition: gtvarlimit.f90:304
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:

◆ gtvarlimit_iiii()

subroutine gtvarlimit_iiii ( type(gt_variable), intent(inout)  var,
integer, intent(in)  dimord,
integer, intent(in), optional  start,
integer, intent(in), optional  count,
integer, intent(in), optional  stride,
logical, intent(out), optional  err 
)

Definition at line 14 of file gtvarlimit.f90.

References dc_trace::beginsub(), dc_trace::dbgmessage(), dc_error::dc_noerr, dc_trace::endsub(), gtdata_internal_map::map_lookup(), gtdata_internal_map::map_set(), and dc_error::storeerror().

14  !
15  !== 入出力範囲の拘束 (数値で指定)
16  !
17  ! 変数 *var* 次元の入出力範囲を拘束します。
18  ! Limit を呼び出した後では Slice でその範囲の外に入出力範囲を
19  ! 設定できなくなります。これにより、変数全体ではなく一部を
20  ! Slice_Next サブルーチンを用いて走査できるようになります。
21  !
22  ! 指定方法は、変数 *var* の *dimord* 番目の次元を基点 *start*,
23  ! 格子総数 *count* 間隔 *stride* に限定します。
24  ! このあとでは、スライス指定の第1格子は *start* 番目の格子を
25  ! 指示することになり、スライス指定での格子数は *count* 個を
26  ! 越えることができなくなり、スライス指定で間引きなしを指定すると
27  ! *stride* 個ごとの指定を指示することになります。
28  !
29  ! エラーが生じた場合、メッセージを出力
30  ! してプログラムは強制終了します。*err* を与えてある場合には
31  ! の引数に .true. が返り、プログラムは終了しません。
32  !
33  ! *Limit* は 2 つのサブルーチンの総称名であり、
34  ! 他にも {gtool4 netCDF 規約}[link:../xref.htm#label-6] の
35  ! 「5.4 コンマ記法」を用いて指定することも可能です。
36  ! 下記のサブルーチンを参照してください。
37  !
38  use gtdata_types, only: gt_variable
40  use dc_error, only: dc_noerr, nf90_einval, storeerror
41  use dc_trace, only: beginsub, endsub, dbgmessage
42  implicit none
43  type(gt_variable), intent(inout):: var
44  integer, intent(in) :: dimord
45  integer, intent(in) , optional :: start, count, stride
46  logical, intent(out), optional :: err
47  type(gt_dimmap), allocatable:: map(:)
48  integer:: iolo, iohi, uilo, uihi, lowerlim, upperlim, dimlo, dimhi
49  integer:: ndims, stat
50 
51  stat = nf90_einval
52  call beginsub('GTVarLimit_iiii', &
53  & 'var%d-dim%d start=%d count=%d stride=%d', &
54  & i=(/var%mapid, dimord, start, count, stride/))
55  ! エラーチェック
56  if (dimord < 1) then
57  print *, "dimord =", dimord, " < 1"
58  goto 999
59  endif
60  if (stride == 0) then
61  print *, "stride == 0"
62  goto 999
63  endif
64  call map_lookup(var, ndims=ndims)
65  if (ndims <= 0) then
66  print *, "ndims =", ndims, " <= 0"
67  goto 999
68  endif
69  if (dimord > ndims) then
70  print *, "dimrod =", dimord, " > ndims =", ndims
71  goto 999
72  endif
73  if (allocated(map)) then
74  deallocate(map)
75  end if
76  allocate(map(ndims))
77  call map_lookup(var, map=map)
78  ! (/lowerlim, upperlim/) は内部格子の範囲 (降順可)
79  lowerlim = min(start, start + (count - 1) * stride)
80  upperlim = max(start, start + (count - 1) * stride)
81  call dimrange(var, dimord, dimlo, dimhi)
82  if (lowerlim < dimlo) then
83  print *, "lowerlim = ", lowerlim, " < dimlo =", dimlo
84  goto 999
85  endif
86  if (upperlim > dimhi) then
87  print *, "upperlim = ", upperlim, " < dimhi =", dimhi
88  goto 999
89  endif
90 
91  call dbgmessage('@ lowerlim=%d upperlim=%d', i=(/lowerlim, upperlim/))
92 
93  ! 入出力範囲を内部格子番号に変えておく
94  uilo = map(dimord)%start
95  iolo = 1 + map(dimord)%step * (uilo - 1) + map(dimord)%offset
96  uihi = map(dimord)%start + (map(dimord)%count - 1) * map(dimord)%stride
97  iohi = 1 + map(dimord)%step * (uihi - 1) + map(dimord)%offset
98 
99  call dbgmessage('@ userindex=%d %d, internal=%d %d', &
100  & i=(/uilo, uihi, iolo, iohi/))
101  call dbgmessage('@ DbgMessage offset %d -> %d step=%d', &
102  & i=(/map(dimord)%offset, (start-1), stride/))
103 
104  ! 制限を課す。offset が変わればユーザ格子番号の意味が変わる
105  map(dimord)%offset = start - 1
106  map(dimord)%allcount = count
107  map(dimord)%step = stride
108 
109  ! 入出力範囲を内部格子番号からユーザ格子番号に戻す
110  uilo = 1 + (iolo - 1 - map(dimord)%offset) / map(dimord)%step
111  uihi = 1 + (iohi - 1 - map(dimord)%offset) / map(dimord)%step
112  call dbgmessage('@ userindex=%d %d', i=(/uilo, uihi/))
113 
114  ! それぞれは制限 [1 .. allcount] の中になければならない
115  uilo = max(1, min(map(dimord)%allcount, uilo))
116  uihi = max(1, min(map(dimord)%allcount, uihi))
117 
118  call dbgmessage('@ userindex=%d %d orig_stride=%d', &
119  & i=(/uilo, uihi, map(dimord)%stride/))
120 
121  ! 元のストライドの符号は無視し、正に固定する
122  map(dimord)%stride = max(1, abs(map(dimord)%stride))
123  map(dimord)%start = min(uilo, uihi)
124  map(dimord)%count = 1 + abs(uihi - uilo) / map(dimord)%stride
125 
126  call map_set(var, map, stat)
127  if (stat /= 0) call dbgmessage("map_set fail")
128 
129 999 continue
130  call storeerror(stat, 'GTVarLimit_iiii', err)
131  call endsub('GTVarLimit_iiii')
subroutine map_set(var, map, stat)
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, 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, public map_lookup(var, vid, map, ndims)
subroutine, public endsub(name, fmt, i, r, d, L, n, c1, c2, c3, ca)
Definition: dc_trace.f90:446
Here is the call graph for this function:

◆ limit_one()

subroutine gtvarlimit::limit_one ( character(len = *), intent(in)  string)

Definition at line 201 of file gtvarlimit.f90.

References dc_url::gt_equal, and region_spec().

Referenced by gtvarlimit(), and gtvarslicec().

201  use dc_url, only: gt_equal
202  use dc_string, only: strieq, stoi
205  character(len = *), intent(in):: string
206  integer:: equal, dimord
207  integer:: start, count, stride, strhead
208  logical:: myerr
209 
210  if (string == '') return
211 
212  strhead = 4
213  if (len(string) < 4) strhead = len(string)
214 
215  if (strieq(string(1:strhead), "IGN:")) then
216  ! 隠蔽型指定子 ign:<dim> または ign:<dim>=<start>
217  equal = index(string, gt_equal)
218  if (equal == 0) then
219  start = 1
220  else
221  start = stoi(string(equal+1: ), default=1)
222  endif
223  dimord = dimname_to_dimord(var, string(5: equal-1))
224  call limit(var, dimord, start, 1, 1, err)
225  call del_dim(var, dimord, myerr)
226  return
227  endif
228 
229  ! 限定型指定子 <dim>=<start>:<finish>:<stride>
230  ! いまは実装がバグっていて <start>:<count>:<stride> になってる
231  !
232  equal = index(string, gt_equal)
233  if (equal == 0) return
234  dimord = dimname_to_dimord(var, string(1: equal-1))
235  if (dimord <= 0) return
236  !
237  call region_spec(dimord, string(equal+1: ), start, count, stride)
238  call limit(var, dimord, start, count, stride, err)
character, parameter, public gt_equal
Definition: dc_url.f90:87
文字型変数の操作.
Definition: dc_string.f90:24
subroutine region_spec(dimord, string, start, count, stride)
Definition: gtvarlimit.f90:245
integer, parameter, public string
文字列を保持する 文字型変数の種別型パラメタ
Definition: dc_types.f90:118
Here is the call graph for this function:
Here is the caller graph for this function:

◆ region_spec()

subroutine gtvarlimit::region_spec ( integer, intent(in)  dimord,
character(len = *), intent(in)  string,
integer, intent(out)  start,
integer, intent(out)  count,
integer, intent(out)  stride 
)

Definition at line 245 of file gtvarlimit.f90.

References dc_url::gt_circumflex, dc_url::gt_colon, dc_types::token, and value_to_index().

Referenced by limit_one().

245  use dc_types, only: token
246  use dc_string, only: index_ofs, stoi
247  use dc_url, only: gt_circumflex, gt_colon
248  use gtdata_internal_map, only: dimrange
249  integer, intent(in):: dimord
250  integer, intent(out):: start, count, stride
251  character(len = *), intent(in):: string
252  integer:: colon, prev_colon, finish, dimlo, dimhi
253  character(len = token):: val(3)
254  continue
255  colon = index(string, gt_colon)
256  if (colon == 0) then
257  ! コロンがない場合は上下端に同じ値
258  val(1) = string(1: )
259  val(2) = val(1)
260  val(3) = ""
261  else
262  val(1) = string(1: colon - 1)
263  prev_colon = colon
264  colon = index_ofs(string, colon + 1, gt_colon)
265  if (colon > 0) then
266  val(2) = string(prev_colon + 1: colon - 1)
267  val(3) = string(colon + 1: )
268  else
269  val(2) = string(prev_colon + 1: )
270  val(3) = ""
271  endif
272  endif
273  if (val(3) == "") val(3) = "^1"
274 
275  if (val(1)(1:1) == gt_circumflex) then
276  start = stoi(val(1)(2: ))
277  else if (val(1) == val(2)) then
278  start = nint(value_to_index(dimord, val(1)))
279  else
280  start = floor(value_to_index(dimord, val(1)))
281  endif
282  if (val(2) == val(1)) then
283  finish = start
284  else if (val(2)(1:1) == gt_circumflex) then
285  finish = stoi(val(2)(2: ))
286  else
287  finish = ceiling(value_to_index(dimord, val(2)))
288  endif
289 
290  call dimrange(var, dimord, dimlo, dimhi)
291  start = min(max(dimlo, start), dimhi)
292  finish = min(max(dimlo, finish), dimhi)
293  count = abs(finish - start) + 1
294 
295  if (val(3)(1:1) == gt_circumflex) then
296  stride = stoi(val(3)(2: ))
297  else
298  stride = stoi(val(3))
299  endif
300  stride = sign(stride, finish - start)
integer, parameter, public token
単語やキーワードを保持する文字型変数の種別型パラメタ
Definition: dc_types.f90:109
character, parameter, public gt_circumflex
Definition: dc_url.f90:89
character, parameter, public gt_colon
Definition: dc_url.f90:83
文字型変数の操作.
Definition: dc_string.f90:24
種別型パラメタを提供します。
Definition: dc_types.f90:49
real function value_to_index(dimord, value)
Definition: gtvarlimit.f90:304
Here is the call graph for this function:
Here is the caller graph for this function:

◆ value_to_index()

real function gtvarlimit::value_to_index ( integer, intent(in)  dimord,
character(len = *), intent(in)  value 
)

Definition at line 304 of file gtvarlimit.f90.

References dc_trace::beginsub(), dc_trace::dbgmessage(), and dc_trace::endsub().

Referenced by region_spec().

304  !
305  ! GTVarLimit の引数 *var* に格納される変数の次元 *dimord*
306  ! に格納されるデータのうち, *value* が格納される
307  ! 格子番号を整数値にして返します.
308  !
309  ! 例えば次元に以下のデータが格納されているとします.
310  !
311  ! 0.05 0.1 0.15 0.20 0.25 0.30
312  !
313  ! この場合, *value* に 0.15 が与えられれば戻り値は 3. となります.
314  ! また *value* に 0.225 が与えられれば戻り値は 4.5 となります.
315  !
316  !
317  use gtdata_types, only: gt_variable
318  use gtdata_generic, only: get, open, close
319  use dc_string, only: stod
320  use dc_trace, only: beginsub, endsub, dbgmessage
321  integer, intent(in):: dimord
322  character(len = *), intent(in):: value
323  type(gt_variable):: axisvar
324  real, pointer:: axisval(:) => null()
325  real:: val
326  integer:: i
327  continue
328 
329  call beginsub('value_to_index', 'var=%d dimord=%d value=%c', &
330  & i=(/var%mapid, dimord/), c1=trim(value))
331 
332  call open(axisvar, var, dimord, count_compact=.true.)
333  nullify(axisval)
334  call get(axisvar, axisval)
335  call close(axisvar)
336  if (.not. associated(axisval)) then
337  result = -1.0
338  return
339  else if (size(axisval) < 2) then
340  result = 1.0
341  goto 900
342  endif
343 
344  val = stod(value)
345 
346  ! call DbgMessage('value=%f axis=(/%*r/)', r=(/val, axisval(:)/), &
347  ! & n=(/size(axisval)/))
348 
349  do, i = 1, size(axisval) - 1
350  if (axisval(i + 1) == axisval(i)) then
351  result = real(i) + 0.5
352  goto 900
353  endif
354  result = i + (val - axisval(i)) / (axisval(i + 1) - axisval(i))
355  if (result <= (i + 1)) goto 900
356  enddo
357 
358 900 continue
359  call endsub('value_to_index', '(%c) = %r', &
360  & c1=trim(value), r=(/result/))
361  deallocate(axisval)
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
文字型変数の操作.
Definition: dc_string.f90:24
subroutine, public endsub(name, fmt, i, r, d, L, n, c1, c2, c3, ca)
Definition: dc_trace.f90:446
Here is the call graph for this function:
Here is the caller graph for this function: