gtvarslice.f90
Go to the documentation of this file.
1 !
2 != 入出力範囲の指定
3 !
4 ! Authors:: Eizi TOYODA, Yasuhiro MORIKAWA
5 ! Version:: $Id: gtvarslice.f90,v 1.3 2009-05-25 09:55:57 morikawa Exp $
6 ! Tag Name:: $Name: $
7 ! Copyright:: Copyright (C) GFD Dennou Club, 2000-2005. All rights reserved.
8 ! License:: See COPYRIGHT[link:../../COPYRIGHT]
9 !
10 ! 以下のサブルーチン、関数は gtdata_generic から gtdata_generic#Slice
11 ! として提供されます。
12 
13 subroutine gtvarslice(var, dimord, start, count, stride)
14  !
15  !== 入出力範囲を数値で指定
16  !
17  ! 変数 *var* の入出力範囲を指定します。
18  !
19  ! 変数 *var* の *dimord* 番目の次元の入出力範囲を *start* から
20  ! *stride* 個おきに *count* 個とします。*start*, *count*,
21  ! *stride* のいずれを省略しても <b>1</b> が仮定されます。成功し
22  ! たか否かを返す引数はありません。仮に指定できない範囲が指定さ
23  ! た場合には、指定範囲を含むなるべく広い範囲を設定します。
24  !
25  ! *Slice* は複数のサブルーチンの総称名であり、
26  ! 他にも文字列や番号で指定する方法があります。
27  !
28  use gtdata_types, only: gt_variable
31  use dc_error, only: nf90_enotvar, storeerror
32  use dc_trace, only: beginsub, endsub, dbgmessage
33  implicit none
34  type(gt_variable), intent(in):: var
35  integer, intent(in):: dimord
36  integer, intent(in), optional:: start
37  integer, intent(in), optional:: count
38  integer, intent(in), optional:: stride
39  type(gt_dimmap), allocatable:: map(:)
40  integer:: vid, maxindex, maxcount, nd, stat
41  logical:: growable_dimension
42 continue
43  call beginsub('GTVarSlice', 'var%%mapid=%d dimord=%d', &
44  & i=(/var%mapid, dimord/))
45  call gtvar_dump(var)
46  call map_lookup(var, vid=vid, ndims=nd)
47  if (vid < 0) then
48  call storeerror(nf90_enotvar, "GTVarSlice")
49  endif
50 
51  if (vid > 0) then
52  call query_growable(vid, growable_dimension)
53  else
54  growable_dimension = .false.
55  endif
56 
57  if (nd == 0) goto 999
58  allocate(map(nd))
59  call map_lookup(var, map=map)
60 
61  if (dimord <= 0 .or. dimord > size(map)) goto 998
62 
63  call dbgmessage('map(dimord): originally start=%d count=%d stride=%d', &
64  & i=(/map(dimord)%start, map(dimord)%count, map(dimord)%stride/))
65  if (.not. growable_dimension) then
66  maxindex = map(dimord)%allcount
67  call dbgmessage('maxindex=%d', i=(/maxindex/))
68  endif
69 
70  if (present(start)) then
71  if (start < 0) then
72  map(dimord)%start = max(1, maxindex + 1 + start)
73  else if (growable_dimension) then
74  map(dimord)%start = max(1, start)
75  else
76  map(dimord)%start = min(maxindex, max(1, start))
77  endif
78  call dbgmessage('start=%d (%d specified)', i=(/map(dimord)%start, start/))
79  endif
80 
81  if (present(stride)) then
82  map(dimord)%stride = stride
83  if (stride == 0) map(dimord)%stride = 1
84  call dbgmessage('stride=%d (%d specified)', &
85  & i=(/map(dimord)%stride, stride/))
86  endif
87 
88  if (present(count)) then
89  map(dimord)%count = abs(count)
90  if (count == 0) map(dimord)%count = 1
91  call dbgmessage('count=%d (%d specified)', &
92  & i=(/map(dimord)%count, count/))
93  endif
94 
95  if (.not. growable_dimension) then
96  maxcount = 1 + (maxindex - map(dimord)%start) / map(dimord)%stride
97  map(dimord)%count = max(1, min(maxcount, map(dimord)%count))
98  call dbgmessage('count=%d ', i=(/map(dimord)%count/))
99  endif
100  call map_set(var, map, stat)
101  if (stat /= 0) goto 998
102 
103  call endsub('GTVarSlice')
104  deallocate(map)
105  return
106 
107 998 continue
108  deallocate(map)
109 999 continue
110  call endsub('GTVarSlice', 'err skipped')
111 end subroutine gtvarslice
112 
113 subroutine gtvarslicec(var, string, err)
114  !
115  !== 入出力範囲を文字列で指定
116  !
117  ! 変数 *var* の入出力範囲を、*string* に応じて指定します。
118  ! *string* には {gtool4 netCDF 規約}[link:../xref.htm#label-6] の
119  ! 「5.4 コンマ記法」に述べられる範囲指定表現を用います。
120  ! 凡例を以下に挙げます。
121  !
122  ! <dim>=<lower>
123  !
124  ! <dim>=<lower>:<upper>
125  !
126  ! <dim>=<lower>:<upper>:<stride>
127  !
128  ! ここで、<tt><dim></tt> は次元番号または次元名であり、
129  ! <tt><lower></tt>, <tt><upper></tt>
130  ! は座標値または "<tt>^</tt>" を前置した格子番号です。
131  ! <tt><stride></tt> は格子数です。
132  !
133  ! 現在 *err* は必ず <tt>.false.</tt> を返すことになっています。
134  !
135  ! *Slice* は複数のサブルーチンの総称名であり、
136  ! 他にも文字列や番号で指定する方法があります。
137  !
138  !
139  !
140  !--
141  ! 変数 var に string によるマップ操作を行う。
142  ! string はコンマで区切られた変換指定の列である。
143  ! 変換指定は領域設定と次元順序変換のどちらかである。
144  ! 領域設定は英数字で始まるもので、<dim>=<lower>, <dim>=<lower>:<upper>,
145  ! <dim>=<lower>:<upper>:<stride> のような形式である。
146  ! ここで、dim は次元番号または次元名であり、<lower>, <upper>
147  ! は ^ を前置した座標即値または格子番号である。
148  ! <stride> は格子数である。
149  ! (未実装) 次元順序変換は = で始まるもので、
150  ! IGN:<dim>=<pos>
151  ! の形態をとる。
152  !++
153  !
154  use gtdata_types, only: gt_variable
155  use gtdata_generic, only: slice
156  use dc_trace, only: beginsub, endsub
157  use dc_url, only: gt_comma
159  type(gt_variable), intent(inout) :: var
160  character(len = *), intent(in) :: string
161  logical, intent(out) :: err
162  integer:: is, ie
163 continue
164  call beginsub('GTVarSliceC', 'var=%d lim=<%c>', &
165  & i=(/var%mapid/), c1=trim(string))
166  call gtvar_dump(var)
167  ! コンマで区切って解釈
168  is = 1
169  do
170  ie = index(string(is: ), gt_comma)
171  if (ie == 0) exit
172  call limit_one(string(is: is+ie-2))
173  is = is + ie
174  if (is > len(string)) exit
175  enddo
176  call limit_one(string(is: ))
177  err = .false.
178  call endsub('GTVarSliceC')
179  return
180 contains
181 
182  subroutine limit_one(string)
183  use dc_url, only: gt_equal
184  use dc_string, only: strieq, stoi
186  character(len = *), intent(in):: string
187  integer:: equal, dimord
188  integer:: start, count, stride
189  logical:: myerr
190 
191  if (string == '') return
192 
193  if (strieq(string(1:4), "IGN:")) then
194  ! 隠蔽型指定子 ign:<dim> または ign:<dim>=<start>
195  equal = index(string, gt_equal)
196  if (equal == 0) then
197  start = 1
198  else
199  start = stoi(string(equal+1: ), default=1)
200  endif
201  dimord = dimname_to_dimord(var, string(5: equal-1))
202  call slice(var, dimord, start, 1, 1)
203  call del_dim(var, dimord, myerr)
204  return
205  endif
206 
207  ! 限定型指定子 <dim>=<start>:<finish>:<stride>
208  !
209  equal = index(string, gt_equal)
210  if (equal == 0) return
211  dimord = dimname_to_dimord(var, string(1: equal-1))
212  if (dimord <= 0) return
213  !
214  call region_spec(dimord, string(equal+1: ), start, count, stride)
215  call slice(var, dimord, start, count, stride)
216  end subroutine limit_one
217 
218  !
219  ! 範囲指定の = のあとを : で区切ってマップにいれる
220  !
221  subroutine region_spec(dimord, string, start, count, stride)
222  use dc_types, only: token
223  use dc_string, only: index_ofs, stoi
224  use dc_url, only: gt_circumflex, gt_colon
225  use gtdata_internal_map, only: dimrange
226  integer, intent(in):: dimord
227  integer, intent(out):: start, count, stride
228  character(len = *), intent(in):: string
229  integer:: colon, prev_colon, finish, dimlo, dimhi
230  character(len = token):: val(3)
231  continue
232  colon = index(string, gt_colon)
233  if (colon == 0) then
234  ! コロンがない場合は上下端に同じ値
235  val(1) = string(1: )
236  val(2) = val(1)
237  val(3) = ""
238  else
239  val(1) = string(1: colon - 1)
240  prev_colon = colon
241  colon = index_ofs(string, colon + 1, gt_colon)
242  if (colon > 0) then
243  val(2) = string(prev_colon + 1: colon - 1)
244  val(3) = string(colon + 1: )
245  else
246  val(2) = string(prev_colon + 1: )
247  val(3) = ""
248  endif
249  endif
250  if (val(3) == "") val(3) = "^1"
251 
252  if (val(1)(1:1) == gt_circumflex) then
253  start = stoi(val(1)(2: ))
254  else if (val(1) == val(2)) then
255  start = nint(value_to_index(dimord, val(1)))
256  else
257  start = floor(value_to_index(dimord, val(1)))
258  endif
259  if (val(2) == val(1)) then
260  finish = start
261  else if (val(2)(1:1) == gt_circumflex) then
262  finish = stoi(val(2)(2: ))
263  else
264  finish = ceiling(value_to_index(dimord, val(2)))
265  endif
266 
267  call dimrange(var, dimord, dimlo, dimhi)
268  start = min(max(dimlo, start), dimhi)
269  finish = min(max(dimlo, finish), dimhi)
270  count = abs(finish - start) + 1
271 
272  if (val(3)(1:1) == gt_circumflex) then
273  stride = stoi(val(3)(2: ))
274  else
275  stride = stoi(val(3))
276  endif
277  stride = sign(stride, finish - start)
278  end subroutine region_spec
279 
280  real function value_to_index(dimord, value) result(result)
281  !
282  ! GTVarSlice の引数 *var* に格納される変数の次元 *dimord*
283  ! に格納されるデータのうち, *value* が格納される
284  ! 格子番号を整数値にして返します.
285  !
286  ! 例えば次元に以下のデータが格納されているとします.
287  !
288  ! 0.05 0.1 0.15 0.20 0.25 0.30
289  !
290  ! この場合, *value* に 0.15 が与えられれば戻り値は 3. となります.
291  ! また *value* に 0.225 が与えられれば戻り値は 4.5 となります.
292  !
293  !
294  use gtdata_types, only: gt_variable
295  use gtdata_generic, only: get, open, close
296  use dc_string, only: stod
297  use dc_trace, only: beginsub, endsub, dbgmessage
298  integer, intent(in):: dimord
299  character(len = *), intent(in):: value
300  type(gt_variable):: axisvar
301  real, pointer:: axisval(:)
302  real:: val
303  integer:: i
304  continue
305  call beginsub('value_to_index', 'var=%d dimord=%d value=%c', &
306  & i=(/var%mapid, dimord/), c1=trim(value))
307 
308  call open(axisvar, var, dimord, count_compact=.true.)
309  call get(axisvar, axisval)
310  call close(axisvar)
311  if (.not. associated(axisval)) then
312  result = -1.0
313  return
314  else if (size(axisval) < 2) then
315  result = 1.0
316  goto 900
317  endif
318 
319  val = stod(value)
320 
321  ! call DbgMessage('value=%f axis=(/%*r/)', r=(/val, axisval(:)/), &
322  ! & n=(/size(axisval)/))
323 
324  do, i = 1, size(axisval) - 1
325  if (axisval(i + 1) == axisval(i)) then
326  result = real(i) + 0.5
327  goto 900
328  endif
329  result = i + (val - axisval(i)) / (axisval(i + 1) - axisval(i))
330  if (result <= (i + 1)) goto 900
331  enddo
332 
333 900 continue
334  call endsub('value_to_index', 'value(%c) =~ index(%r)', &
335  & c1=trim(value), r=(/result/))
336  deallocate(axisval)
337  end function value_to_index
338 
339 end subroutine gtvarslicec
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 gtvarslicec(var, string, err)
Definition: gtvarslice.f90:114
subroutine map_set(var, map, stat)
subroutine, public storeerror(number, where, err, cause_c, cause_i)
Definition: dc_error.f90:830
subroutine gtvarslice(var, dimord, start, count, stride)
Definition: gtvarslice.f90:14
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
subroutine, public map_lookup(var, vid, map, ndims)
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
subroutine, public query_growable(vid, result)