gtdata_internal_map.f90
Go to the documentation of this file.
1 != gtool 変数表に関する各種手続き
2 !
3 ! Authors:: Eizi TOYODA, Yasuhiro MORIKAWA
4 ! Version:: $Id: gtdata_internal_map.f90,v 1.3 2010-04-11 14:13:50 morikawa Exp $
5 ! Tag Name:: $Name: $
6 ! Copyright:: Copyright (C) GFD Dennou Club, 2000-2005. All rights reserved.
7 ! License:: See COPYRIGHT[link:../../COPYRIGHT]
8 !
9 
11  !
12  !== gtool 変数表
13  !
14  ! gtool 変数というのは実はマップ表のキーとなる整数ハンドルである。
15  ! マップ表 maptab には実体表のエントリ番号と次元書き換え/イテレータ
16  ! の表が載っている。
17  ! このレベルにおける参照カウントは作らないことにする。つまり、
18  ! マップ表と実体表は一対一対応するし、
19  ! ユーザがハンドルをコピーするのは勝手である。
20  ! もちろんユーザには必ずただ1回
21  ! 当該ハンドルを close すなわち maptabdelete する義務がある。
22 
23  use dc_types, only: string
26  ! これらは, 外部から gtdata_internal_vartable を直接呼ばないようにするため,
27  ! gtdata_internal_map を介して公開する.
28 
29  implicit none
30 
31  type gt_dimmap
32  ! 次元書き換え表
33  integer:: dimno
34  ! 正ならば実体変数の次元番号, 他変数参照時は非正値
35  character(len=STRING):: url
36  ! 次元変数の url
37  integer:: offset
38  !=== 実体と gtool ユーザの格子番号対応
39  !
40  ! ユーザからみて 1..allcount が実体の
41  ! (1..allcount) * step + offset に写像される。
42  ! これらの値の変更の際は実体変数の許容添字範囲および
43  ! 成長可能性を確認する必要あり。
44  ! start 値に対するオフセット。
45  integer:: step
46  ! 1 か -1 の値をとることが期待される。
47  integer:: allcount
48  ! 見掛けの格子番号上限: start + count * stride <= allcount
49  integer:: start
50  !=== イテレータ本体
51  ! 入出力範囲は (start:start+count*stride:stride) である。
52  ! イテレータ start 値
53  integer:: count
54  ! イテレータ count 値
55  integer:: stride
56  ! イテレータ stride 値
57  logical:: scalar
58  ! スカラー変数である場合, .true.
59  end type gt_dimmap
60 
62  ! マップ表の型
63  integer:: vid
64  integer:: ndims
65  type(gt_dimmap), pointer:: map(:)
66  end type map_table_entry
67 
68  type(map_table_entry), save, target, allocatable:: maptab(:)
69  integer, parameter:: maptab_init_size = 16
70 
74  private:: maptab, maptab_init_size
75 
76  interface dimrange
77  module procedure dimrange_by_dimno
78  end interface
79 
80 contains
81 
82  subroutine dimrange_by_dimno(var, dimno, dimlo, dimhi)
83  ! 変数と次元番号を指定して、当該次元の内部的添字番号範囲を得る
84  use gtdata_types, only: gt_variable
85  use gtdata_generic, only: open, close
87  type(gt_variable), intent(in):: var
88  integer, intent(in):: dimno
89  integer, intent(out):: dimlo, dimhi
90  type(gt_variable):: dimvar
91  integer:: vid
92  call open(dimvar, var, dimno, count_compact=.true.)
93  call map_lookup(dimvar, vid=vid)
94  call dimrange(vid, dimlo, dimhi)
95  call close(dimvar)
96  end subroutine dimrange_by_dimno
97 
98  subroutine map_dup(var, source_var)
99  ! 変数 source_var の複写 var を作成する
100  use gtdata_types, only: gt_variable
102  use dc_trace, only: dbgmessage
103  type(gt_variable), intent(out):: var
104  type(gt_variable), intent(in):: source_var
105  integer:: vid, mid1, mid2, vid2, nd, class, cid
106  call map_lookup(source_var, vid=vid)
107  if (vid < 0) then
108  var = gt_variable(-1)
109  return
110  endif
111  if (vid == 0) then
112  vid2 = 0
113  else
114  call vartablelookup(vid, class=class, cid=cid)
115  call vartableadd(vid2, class, cid)
116  endif
117  call maptabadd(var%mapid, vid2)
118  mid1 = source_var%mapid
119  mid2 = var%mapid
120  maptab(mid2)%ndims = maptab(mid1)%ndims
121  if (associated(maptab(mid1)%map)) then
122  nd = size(maptab(mid1)%map)
123  allocate(maptab(mid2)%map(nd))
124  maptab(mid2)%map(1:nd) = maptab(mid1)%map(1:nd)
125  else
126  nullify(maptab(mid2)%map)
127  endif
128  call dbgmessage('map_dup mapid(%d from %d) vid(%d from %d)', &
129  & i=(/mid2, mid1, maptab(mid2)%vid, maptab(mid1)%vid/))
130  end subroutine map_dup
131 
132  subroutine map_create(var, class, cid, ndims, allcount, stat)
133  ! 変数 var を作成する。内部種別 class, 内部識別子 cid,
134  ! 外見的次元数 ndims, 外見的次元長 allcount(:) を与える。
135  ! オフセットゼロを仮定して諸元の初期化が行われる。
136  use gtdata_types, only: gt_variable
138  use dc_error, only: nf90_enotvar, gt_enomoredims, dc_noerr
139  type(gt_variable), intent(out):: var
140  integer, intent(in):: class, cid, ndims, allcount(:)
141  integer, intent(out):: stat
142  type(gt_dimmap), pointer:: map(:)
143  integer:: vid, i
144  continue
145 
146  stat = dc_noerr
147  if ( ndims < 0 ) then
148  stat = gt_enomoredims
149  goto 999
150  end if
151  call vartableadd(vid, class, cid)
152  call maptabadd(var%mapid, vid)
153  if (ndims > 0) then
154  call map_allocate(map, ndims)
155  maptab(var%mapid)%ndims = ndims
156  maptab(var%mapid)%map => map
157 
158  do, i = 1, ndims
159  map(i)%dimno = i
160  map(i)%allcount = allcount(i)
161  map(i)%count = allcount(i)
162  map(i)%offset = 0
163  map(i)%start = 1
164  map(i)%step = 1
165  map(i)%stride = 1
166  map(i)%scalar = .false.
167  enddo
168  else
169  ! スカラー変数 (ndims = 0) の場合
170  call map_allocate(map, 1)
171  maptab(var%mapid)%ndims = 0
172  maptab(var%mapid)%map => map
173  map(1)%dimno = 1
174  map(1)%allcount = 1
175  map(1)%count = 1
176  map(1)%offset = 0
177  map(1)%start = 1
178  map(1)%step = 1
179  map(1)%stride = 1
180  map(1)%scalar = .true.
181  end if
182 
183 999 continue
184  return
185  end subroutine map_create
186 
187  subroutine maptabadd(mapid, vid)
188  ! すでに実体表に追加されたエントリ番号 vid を指定して、
189  ! マップ表にエントリを追加する。
190  integer, intent(out):: mapid
191  integer, intent(in):: vid
192  type(map_table_entry), allocatable:: tmp_maptab(:)
193  integer:: i, n
194  ! 必要なら初期確保
195  if (.not. allocated(maptab)) then
196  allocate(maptab(maptab_init_size))
197  maptab(:)%vid = vid_invalid
198  do, n = 1, maptab_init_size
199  nullify(maptab(n)%map)
200  enddo
201  endif
202  ! 空き地があればそこに割り当て
203  do, i = 1, size(maptab)
204  if (maptab(i)%vid == vid_invalid) then
205  mapid = i
206  maptab(mapid)%vid = vid
207  return
208  endif
209  enddo
210  ! 空き地はなかったのだから倍幅確保
211  n = size(maptab)
212  allocate(tmp_maptab(n))
213  tmp_maptab(:) = maptab(:)
214  deallocate(maptab)
215  allocate(maptab(n * 2))
216  ! 確保したところはクリア
217  maptab(1:n) = tmp_maptab(1:n)
218  do, i = n + 1, (2 * size(tmp_maptab))
219  maptab(i)%vid = vid_invalid
220  nullify(maptab(i)%map)
221  enddo
222  deallocate(tmp_maptab)
223  mapid = n + 1
224  maptab(mapid)%vid = vid
225  end subroutine maptabadd
226 
227  subroutine maptabdelete(var, err)
228  ! 変数 var をマップ表から削除する。
229  ! 実体表には手をつけない。
230  use dc_error, only: nf90_enotvar, storeerror, dc_noerr
231  use gtdata_types, only: gt_variable
232  use dc_trace, only: dbgmessage
233  implicit none
234  type(gt_variable), intent(in):: var
235  logical, intent(out), optional:: err
236  integer:: mapid
237  mapid = var%mapid
238  if (.not. allocated(maptab)) goto 999
239  if (mapid <= 0 .or. mapid > size(maptab)) goto 999
240  if (maptab(mapid)%vid == vid_invalid) goto 999
241  maptab(mapid)%vid = vid_invalid
242  if (associated(maptab(mapid)%map)) deallocate(maptab(mapid)%map)
243  call storeerror(dc_noerr, 'maptabdelete', err)
244  call dbgmessage('gtdata_internal_map table %d deleted', i=(/mapid/))
245  return
246 999 continue
247  call storeerror(nf90_enotvar, 'maptabdelete', err)
248  end subroutine maptabdelete
249 
250  subroutine map_lookup(var, vid, map, ndims)
251  ! 同じファイル番号の変数表の中身を返す
252  use gtdata_types, only: gt_variable
253  type(gt_variable), intent(in):: var
254  integer, intent(out), optional:: vid
255  type(gt_dimmap), intent(out), optional:: map(:)
256  integer, intent(out), optional:: ndims
257  if (.not. allocated(maptab)) goto 999
258  if (var%mapid <= 0 .or. var%mapid > size(maptab)) goto 999
259  if (maptab(var%mapid)%vid == vid_invalid) goto 999
260  if (present(vid)) vid = maptab(var%mapid)%vid
261  if (present(map)) map(:) = maptab(var%mapid)%map(1:size(map))
262  if (present(ndims)) ndims = maptab(var%mapid)%ndims
263  return
264 999 continue
265  if (present(vid)) vid = vid_invalid
266  if (present(map)) then
267  map(:)%dimno = -1
268  map(:)%url = " "
269  endif
270  if (present(ndims)) ndims = 0
271  end subroutine map_lookup
272 
273  subroutine map_set(var, map, stat)
274  ! 同じファイル番号の変数表の値を設定する
275  use gtdata_types, only: gt_variable
276  use dc_error, only: nf90_enotvar, gt_enomoredims, dc_noerr
277  type(gt_variable), intent(in):: var
278  type(gt_dimmap), intent(in):: map(:)
279  integer, intent(out):: stat
280  if (.not. allocated(maptab)) goto 999
281  if (var%mapid <= 0 .or. var%mapid > size(maptab)) goto 999
282  if (maptab(var%mapid)%vid == vid_invalid) goto 999
283  if (size(map) > size(maptab(var%mapid)%map)) then
284  stat = gt_enomoredims
285  return
286  endif
287  maptab(var%mapid)%map(1:size(map)) = map(:)
288  stat = dc_noerr
289  return
290 999 continue
291  stat = nf90_enotvar
292  end subroutine map_set
293 
294 
295  subroutine var_class(var, class, cid)
296  ! 変数 var を指定して、内部種別 class, 内部識別子 cid を得る。
297  use gtdata_types, only: gt_variable
299  type(gt_variable), intent(in):: var
300  integer, intent(out), optional:: class, cid
301  integer:: vid
302  call map_lookup(var, vid=vid)
303  call vartablelookup(vid, class=class, cid=cid)
304  end subroutine var_class
305 
306  subroutine map_set_ndims(var, ndims, stat)
307  ! 変数 var の次元数を ndims に変える。
308  use gtdata_types, only: gt_variable
310  use dc_error, only: nf90_enotvar, gt_enomoredims, dc_noerr
311  type(gt_variable), intent(in):: var
312  integer, intent(in):: ndims
313  integer, intent(out):: stat
314  integer:: vid
315  call map_lookup(var, vid=vid)
316  if (vid == vid_invalid) then
317  stat = nf90_enotvar
318  return
319  endif
320  if (.not. associated(maptab(var%mapid)%map)) then
321  if (ndims == 0) then
322  stat = dc_noerr
323  maptab(var%mapid)%ndims = 0
324  else
325  stat = gt_enomoredims
326  endif
327  else
328  if (ndims > size(maptab(var%mapid)%map)) then
329  stat = gt_enomoredims
330  else
331  stat = dc_noerr
332  maptab(var%mapid)%ndims = ndims
333  endif
334  endif
335  end subroutine map_set_ndims
336 
337  subroutine map_set_rank(var, rank, stat)
338  ! 変数 var のランク(非縮退次元数)を rank に減らすように
339  ! count 値を1に減らす。ランクを増やすことや外見次元数の操作はしない。
340  use gtdata_types, only: gt_variable
342  use dc_error, only: nf90_enotvar, gt_enomoredims, dc_noerr
343  type(gt_variable), intent(in):: var
344  integer, intent(in):: rank
345  integer, intent(out):: stat
346  type(gt_dimmap), pointer:: tmpmap(:)
347  integer:: ndims
348  integer:: vid, nd
349  call map_lookup(var, vid, ndims=ndims)
350  if (vid == vid_invalid) then
351  stat = nf90_enotvar
352  return
353  endif
354  if (ndims < rank) then
355  stat = gt_enomoredims
356  return
357  endif
358  tmpmap => maptab(var%mapid)%map
359  do, nd = ndims, 1, -1
360  if (count(tmpmap(1:ndims)%count > 1) <= rank) exit
361  tmpmap(nd)%count = 1
362  enddo
363  stat = dc_noerr
364  end subroutine map_set_rank
365 
366  subroutine map_to_internal_specs(var, specs, ndims)
367  ! マップ表から netCDF の引数にふさわしい start, count, stride, imap
368  ! を作成する。ただし、stride が負になるばあいは対策されていない。
369  ! (暫定的に gdncvarget/gdncvarput が対応している)
370  use gtdata_types, only: gt_variable
371  use gtdata_internal_vartable, only: num_dimensions => ndims
372  type(gt_variable), intent(in):: var
373  integer, pointer:: specs(:, :)
374  integer, intent(out), optional:: ndims
375  type(gt_dimmap), pointer:: it
376  integer:: vid, i, j, imap, internal_ndims
377  integer:: external_ndims
378  continue
379  call map_lookup(var, vid, ndims=external_ndims)
380  internal_ndims = num_dimensions(vid)
381  if (present(ndims)) ndims = internal_ndims
382  allocate(specs(max(1, internal_ndims), 4))
383  specs(:, 1) = 1
384  specs(:, 2) = 1
385  specs(:, 3) = 1
386  specs(:, 4) = 0
387  imap = 1
388  do, i = 1, size(maptab(var%mapid)%map)
389  it => maptab(var%mapid)%map(i)
390  j = it%dimno
391  if (j > 0 .and. j <= internal_ndims) then
392  specs(j, 1) = it%start + it%offset
393  specs(j, 2) = it%count
394  if (i > external_ndims) specs(j, 2) = 1
395  specs(j, 3) = it%stride * it%step
396  specs(j, 4) = imap
397  endif
398  imap = imap * it%count
399  enddo
400  end subroutine map_to_internal_specs
401 
402  subroutine map_allocate(map, ndims)
403  ! 次元表エントリに ndims 個のエントリを割り付け初期化する。
404  type(gt_dimmap), pointer:: map(:)
405  integer, intent(in):: ndims
406  if (ndims <= 0) then
407  nullify(map)
408  return
409  endif
410  allocate(map(1:ndims))
411  map(1:ndims)%dimno = -1
412  map(1:ndims)%url = ' '
413  map(1:ndims)%allcount = 0
414  map(1:ndims)%offset = 0
415  map(1:ndims)%step = 1
416  map(1:ndims)%start = 1
417  map(1:ndims)%count = 0
418  map(1:ndims)%stride = 1
419  map(1:ndims)%scalar = .false.
420  end subroutine map_allocate
421 
422  subroutine map_apply(var, map)
423  ! 変数 var にマップ表 map を組み合わせる
424  use gtdata_types, only: gt_variable
425  type(gt_variable), intent(inout):: var
426  type(gt_dimmap), pointer:: map(:)
427  type(gt_dimmap), pointer:: tmpmap(:), varmap
428  integer:: i, nd
429  nd = size(map)
430  allocate(tmpmap(nd))
431  do, i = 1, nd
432  tmpmap(i)%allcount = map(i)%allcount
433  tmpmap(i)%count = map(i)%count
434  if (map(i)%dimno > 0) then
435  varmap => maptab(var%mapid)%map(map(i)%dimno)
436  tmpmap(i)%url = varmap%url
437  tmpmap(i)%dimno = varmap%dimno
438  tmpmap(i)%offset = varmap%offset + map(i)%offset
439  tmpmap(i)%step = varmap%step * map(i)%step
440  else
441  tmpmap(i)%url = map(i)%url
442  tmpmap(i)%dimno = 0
443  tmpmap(i)%offset = map(i)%offset
444  tmpmap(i)%step = map(i)%step
445  endif
446  enddo
447  deallocate(map)
448  map => tmpmap
449  end subroutine map_apply
450 
451  subroutine map_resize(var, ndims)
452  ! 変数 var の次元表の大きさを変える
453  use gtdata_types, only: gt_variable
454  type(gt_variable), intent(in):: var
455  integer, intent(in):: ndims
456  type(gt_dimmap), pointer:: newmap(:)
457  type(gt_dimmap), pointer:: tmpmap(:)
458  integer:: n
459  if (associated(maptab(var%mapid)%map)) then
460  tmpmap => maptab(var%mapid)%map
461  call map_allocate(newmap, ndims)
462  n = min(size(tmpmap), ndims)
463  newmap(1:n) = tmpmap(1:n)
464  deallocate(tmpmap)
465  maptab(var%mapid)%map => newmap
466  newmap(n+1:ndims)%dimno = -1
467  newmap(n+1:ndims)%url = ' '
468  newmap(n+1:ndims)%allcount = 0
469  newmap(n+1:ndims)%offset = 0
470  newmap(n+1:ndims)%step = 1
471  newmap(n+1:ndims)%start = 1
472  newmap(n+1:ndims)%count = 0
473  newmap(n+1:ndims)%stride = 1
474  else
475  call map_allocate(maptab(var%mapid)%map, ndims)
476  n = 1
477  endif
478  end subroutine map_resize
479 
480  subroutine gtvar_dump(var)
481  ! 変数のプロパティを出力
482  use gtdata_types, only: gt_variable
484  use dc_trace, only: debug, dbgmessage
485  type(gt_variable), intent(in):: var
486  integer:: idim, imap
487  logical:: dbg_mode
488  continue
489  call debug( dbg_mode )
490  if (.not. dbg_mode) return
491  imap = var%mapid
492  if (imap < 1 .or. imap > size(maptab)) then
493  call dbgmessage('[gt_variable %d: invalid id]', i=(/imap/))
494  return
495  endif
496  if (associated(maptab(imap)%map)) then
497  call dbgmessage('[gt_variable %d: ndims=%d, map.size=%d]', &
498  & i=(/imap, maptab(imap)%ndims, size(maptab(imap)%map)/))
499  do, idim = 1, size(maptab(imap)%map)
500  call dbgmessage('[dim%d dimno=%d ofs=%d step=%d' &
501  &// ' all=%d start=%d count=%d stride=%d url=%c]', &
502  & c1=trim(maptab(imap)%map(idim)%url), &
503  & i=(/idim, maptab(imap)%map(idim)%dimno, &
504  & maptab(imap)%map(idim)%offset, &
505  & maptab(imap)%map(idim)%step, &
506  & maptab(imap)%map(idim)%allcount, &
507  & maptab(imap)%map(idim)%start, &
508  & maptab(imap)%map(idim)%count, &
509  & maptab(imap)%map(idim)%stride/))
510  enddo
511  else
512  call dbgmessage('[gt_variable %d: ndims=%d, map=null]', &
513  & i=(/imap, maptab(imap)%ndims/))
514  endif
515  call vartable_dump(maptab(imap)%vid)
516  end subroutine gtvar_dump
517 
518  integer function dimord_skip_compact(dimord, map) result(result)
519  ! 次元表の中で非縮退次元だけを数えた次元番号 dimord の次元を
520  ! 特定し、外部向けの次元番号を返す。
521  use dc_trace, only: dbgmessage
522  integer, intent(in):: dimord
523  type(gt_dimmap), intent(in):: map(:)
524  integer:: nd, id
525  result = -1
526  nd = 0
527  do, id = 1, size(map)
528  if (map(id)%count < 2) cycle
529  nd = nd + 1
530  if (nd < dimord) cycle
531  result = id
532  call dbgmessage('compact dim skip: %d <= %d', i=(/result, dimord/))
533  exit
534  enddo
535  end function dimord_skip_compact
536 
537 end module gtdata_internal_map
subroutine map_apply(var, map)
subroutine, public map_to_internal_specs(var, specs, ndims)
subroutine dimrange_by_dimno(var, dimno, dimlo, dimhi)
subroutine, public vartable_dump(vid)
integer, parameter, public vtb_class_netcdf
integer, parameter, public vid_invalid
subroutine map_dup(var, source_var)
subroutine, public maptabdelete(var, err)
integer function dimord_skip_compact(dimord, map)
subroutine, public map_create(var, class, cid, ndims, allcount, stat)
subroutine map_set(var, map, stat)
subroutine, public storeerror(number, where, err, cause_c, cause_i)
Definition: dc_error.f90:830
subroutine map_set_ndims(var, ndims, stat)
integer, parameter, public dc_noerr
Definition: dc_error.f90:509
integer function, public ndims(vid)
subroutine map_set_rank(var, rank, stat)
integer, parameter, public vtb_class_unused
subroutine, public vartablelookup(vid, class, cid)
Provides kind type parameter values.
Definition: dc_types.f90:49
subroutine, public maptabadd(mapid, vid)
subroutine, public map_lookup(var, vid, map, ndims)
subroutine, public vartableadd(vid, class, cid)
integer, parameter, public gt_enomoredims
Definition: dc_error.f90:528
subroutine map_allocate(map, ndims)
subroutine, public var_class(var, class, cid)
subroutine map_resize(var, ndims)
integer, parameter, public string
Character length for string.
Definition: dc_types.f90:118