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 
73  public:: vid_invalid, &
75  private:: maptab, maptab_init_size
76 
77  interface dimrange
78  module procedure dimrange_by_dimno
79  end interface
80 
81 contains
82 
83  subroutine dimrange_by_dimno(var, dimno, dimlo, dimhi)
84  ! 変数と次元番号を指定して、当該次元の内部的添字番号範囲を得る
85  use gtdata_types, only: gt_variable
86  use gtdata_generic, only: open, close
88  type(gt_variable), intent(in):: var
89  integer, intent(in):: dimno
90  integer, intent(out):: dimlo, dimhi
91  type(gt_variable):: dimvar
92  integer:: vid
93  call open(dimvar, var, dimno, count_compact=.true.)
94  call map_lookup(dimvar, vid=vid)
95  call dimrange(vid, dimlo, dimhi)
96  call close(dimvar)
97  end subroutine dimrange_by_dimno
98 
99  subroutine map_dup(var, source_var)
100  ! 変数 source_var の複写 var を作成する
101  use gtdata_types, only: gt_variable
103  use dc_trace, only: dbgmessage
104  type(gt_variable), intent(out):: var
105  type(gt_variable), intent(in):: source_var
106  integer:: vid, mid1, mid2, vid2, nd, class, cid
107  call map_lookup(source_var, vid=vid)
108  if (vid < 0) then
109  var = gt_variable(-1)
110  return
111  endif
112  if (vid == 0) then
113  vid2 = 0
114  else
115  call vartablelookup(vid, class=class, cid=cid)
116  call vartableadd(vid2, class, cid)
117  endif
118  call maptabadd(var%mapid, vid2)
119  mid1 = source_var%mapid
120  mid2 = var%mapid
121  maptab(mid2)%ndims = maptab(mid1)%ndims
122  if (associated(maptab(mid1)%map)) then
123  nd = size(maptab(mid1)%map)
124  allocate(maptab(mid2)%map(nd))
125  maptab(mid2)%map(1:nd) = maptab(mid1)%map(1:nd)
126  else
127  nullify(maptab(mid2)%map)
128  endif
129  call dbgmessage('map_dup mapid(%d from %d) vid(%d from %d)', &
130  & i=(/mid2, mid1, maptab(mid2)%vid, maptab(mid1)%vid/))
131  end subroutine map_dup
132 
133  subroutine map_create(var, class, cid, ndims, allcount, stat)
134  ! 変数 var を作成する。内部種別 class, 内部識別子 cid,
135  ! 外見的次元数 ndims, 外見的次元長 allcount(:) を与える。
136  ! オフセットゼロを仮定して諸元の初期化が行われる。
137  use gtdata_types, only: gt_variable
139  use dc_error, only: nf90_enotvar, gt_enomoredims, dc_noerr
140  type(gt_variable), intent(out):: var
141  integer, intent(in):: class, cid, ndims, allcount(:)
142  integer, intent(out):: stat
143  type(gt_dimmap), pointer:: map(:)
144  integer:: vid, i
145  continue
146 
147  stat = dc_noerr
148  if ( ndims < 0 ) then
149  stat = gt_enomoredims
150  goto 999
151  end if
152  call vartableadd(vid, class, cid)
153  call maptabadd(var%mapid, vid)
154  if (ndims > 0) then
155  call map_allocate(map, ndims)
156  maptab(var%mapid)%ndims = ndims
157  maptab(var%mapid)%map => map
158 
159  do, i = 1, ndims
160  map(i)%dimno = i
161  map(i)%allcount = allcount(i)
162  map(i)%count = allcount(i)
163  map(i)%offset = 0
164  map(i)%start = 1
165  map(i)%step = 1
166  map(i)%stride = 1
167  map(i)%scalar = .false.
168  enddo
169  else
170  ! スカラー変数 (ndims = 0) の場合
171  call map_allocate(map, 1)
172  maptab(var%mapid)%ndims = 0
173  maptab(var%mapid)%map => map
174  map(1)%dimno = 1
175  map(1)%allcount = 1
176  map(1)%count = 1
177  map(1)%offset = 0
178  map(1)%start = 1
179  map(1)%step = 1
180  map(1)%stride = 1
181  map(1)%scalar = .true.
182  end if
183 
184 999 continue
185  return
186  end subroutine map_create
187 
188  subroutine maptabadd(mapid, vid)
189  ! すでに実体表に追加されたエントリ番号 vid を指定して、
190  ! マップ表にエントリを追加する。
191  integer, intent(out):: mapid
192  integer, intent(in):: vid
193  type(map_table_entry), allocatable:: tmp_maptab(:)
194  integer:: i, n
195  ! 必要なら初期確保
196  if (.not. allocated(maptab)) then
197  allocate(maptab(maptab_init_size))
198  maptab(:)%vid = vid_invalid
199  do, n = 1, maptab_init_size
200  nullify(maptab(n)%map)
201  enddo
202  endif
203  ! 空き地があればそこに割り当て
204  do, i = 1, size(maptab)
205  if (maptab(i)%vid == vid_invalid) then
206  mapid = i
207  maptab(mapid)%vid = vid
208  return
209  endif
210  enddo
211  ! 空き地はなかったのだから倍幅確保
212  n = size(maptab)
213  allocate(tmp_maptab(n))
214  tmp_maptab(:) = maptab(:)
215  deallocate(maptab)
216  allocate(maptab(n * 2))
217  ! 確保したところはクリア
218  maptab(1:n) = tmp_maptab(1:n)
219  do, i = n + 1, (2 * size(tmp_maptab))
220  maptab(i)%vid = vid_invalid
221  nullify(maptab(i)%map)
222  enddo
223  deallocate(tmp_maptab)
224  mapid = n + 1
225  maptab(mapid)%vid = vid
226  end subroutine maptabadd
227 
228  subroutine maptabdelete(var, err)
229  ! 変数 var をマップ表から削除する。
230  ! 実体表には手をつけない。
231  use dc_error, only: nf90_enotvar, storeerror, dc_noerr
232  use gtdata_types, only: gt_variable
233  use dc_trace, only: dbgmessage
234  implicit none
235  type(gt_variable), intent(in):: var
236  logical, intent(out), optional:: err
237  integer:: mapid
238  mapid = var%mapid
239  if (.not. allocated(maptab)) goto 999
240  if (mapid <= 0 .or. mapid > size(maptab)) goto 999
241  if (maptab(mapid)%vid == vid_invalid) goto 999
242  maptab(mapid)%vid = vid_invalid
243  if (associated(maptab(mapid)%map)) deallocate(maptab(mapid)%map)
244  call storeerror(dc_noerr, 'maptabdelete', err)
245  call dbgmessage('gtdata_internal_map table %d deleted', i=(/mapid/))
246  return
247 999 continue
248  call storeerror(nf90_enotvar, 'maptabdelete', err)
249  end subroutine maptabdelete
250 
251  subroutine map_lookup(var, vid, map, ndims)
252  ! 同じファイル番号の変数表の中身を返す
253  use gtdata_types, only: gt_variable
254  type(gt_variable), intent(in):: var
255  integer, intent(out), optional:: vid
256  type(gt_dimmap), intent(out), optional:: map(:)
257  integer, intent(out), optional:: ndims
258  if (.not. allocated(maptab)) goto 999
259  if (var%mapid <= 0 .or. var%mapid > size(maptab)) goto 999
260  if (maptab(var%mapid)%vid == vid_invalid) goto 999
261  if (present(vid)) vid = maptab(var%mapid)%vid
262  if (present(map)) map(:) = maptab(var%mapid)%map(1:size(map))
263  if (present(ndims)) ndims = maptab(var%mapid)%ndims
264  return
265 999 continue
266  if (present(vid)) vid = vid_invalid
267  if (present(map)) then
268  map(:)%dimno = -1
269  map(:)%url = " "
270  endif
271  if (present(ndims)) ndims = 0
272  end subroutine map_lookup
273 
274  subroutine map_set(var, map, stat)
275  ! 同じファイル番号の変数表の値を設定する
276  use gtdata_types, only: gt_variable
277  use dc_error, only: nf90_enotvar, gt_enomoredims, dc_noerr
278  type(gt_variable), intent(in):: var
279  type(gt_dimmap), intent(in):: map(:)
280  integer, intent(out):: stat
281  if (.not. allocated(maptab)) goto 999
282  if (var%mapid <= 0 .or. var%mapid > size(maptab)) goto 999
283  if (maptab(var%mapid)%vid == vid_invalid) goto 999
284  if (size(map) > size(maptab(var%mapid)%map)) then
285  stat = gt_enomoredims
286  return
287  endif
288  maptab(var%mapid)%map(1:size(map)) = map(:)
289  stat = dc_noerr
290  return
291 999 continue
292  stat = nf90_enotvar
293  end subroutine map_set
294 
295 
296  subroutine var_class(var, class, cid)
297  ! 変数 var を指定して、内部種別 class, 内部識別子 cid を得る。
298  use gtdata_types, only: gt_variable
300  type(gt_variable), intent(in):: var
301  integer, intent(out), optional:: class, cid
302  integer:: vid
303  call map_lookup(var, vid=vid)
304  call vartablelookup(vid, class=class, cid=cid)
305  end subroutine var_class
306 
307  subroutine map_set_ndims(var, ndims, stat)
308  ! 変数 var の次元数を ndims に変える。
309  use gtdata_types, only: gt_variable
311  use dc_error, only: nf90_enotvar, gt_enomoredims, dc_noerr
312  type(gt_variable), intent(in):: var
313  integer, intent(in):: ndims
314  integer, intent(out):: stat
315  integer:: vid
316  call map_lookup(var, vid=vid)
317  if (vid == vid_invalid) then
318  stat = nf90_enotvar
319  return
320  endif
321  if (.not. associated(maptab(var%mapid)%map)) then
322  if (ndims == 0) then
323  stat = dc_noerr
324  maptab(var%mapid)%ndims = 0
325  else
326  stat = gt_enomoredims
327  endif
328  else
329  if (ndims > size(maptab(var%mapid)%map)) then
330  stat = gt_enomoredims
331  else
332  stat = dc_noerr
333  maptab(var%mapid)%ndims = ndims
334  endif
335  endif
336  end subroutine map_set_ndims
337 
338  subroutine map_set_rank(var, rank, stat)
339  ! 変数 var のランク(非縮退次元数)を rank に減らすように
340  ! count 値を1に減らす。ランクを増やすことや外見次元数の操作はしない。
341  use gtdata_types, only: gt_variable
343  use dc_error, only: nf90_enotvar, gt_enomoredims, dc_noerr
344  type(gt_variable), intent(in):: var
345  integer, intent(in):: rank
346  integer, intent(out):: stat
347  type(gt_dimmap), pointer:: tmpmap(:)
348  integer:: ndims
349  integer:: vid, nd
350  call map_lookup(var, vid, ndims=ndims)
351  if (vid == vid_invalid) then
352  stat = nf90_enotvar
353  return
354  endif
355  if (ndims < rank) then
356  stat = gt_enomoredims
357  return
358  endif
359  tmpmap => maptab(var%mapid)%map
360  do, nd = ndims, 1, -1
361  if (count(tmpmap(1:ndims)%count > 1) <= rank) exit
362  tmpmap(nd)%count = 1
363  enddo
364  stat = dc_noerr
365  end subroutine map_set_rank
366 
367  subroutine map_to_internal_specs(var, specs, ndims)
368  ! マップ表から netCDF の引数にふさわしい start, count, stride, imap
369  ! を作成する。ただし、stride が負になるばあいは対策されていない。
370  ! (暫定的に gdncvarget/gdncvarput が対応している)
371  use gtdata_types, only: gt_variable
372  use gtdata_internal_vartable, only: num_dimensions => ndims
373  type(gt_variable), intent(in):: var
374  integer, pointer:: specs(:, :)
375  integer, intent(out), optional:: ndims
376  type(gt_dimmap), pointer:: it
377  integer:: vid, i, j, imap, internal_ndims
378  integer:: external_ndims
379  continue
380  call map_lookup(var, vid, ndims=external_ndims)
381  internal_ndims = num_dimensions(vid)
382  if (present(ndims)) ndims = internal_ndims
383  allocate(specs(max(1, internal_ndims), 4))
384  specs(:, 1) = 1
385  specs(:, 2) = 1
386  specs(:, 3) = 1
387  specs(:, 4) = 0
388  imap = 1
389  do, i = 1, size(maptab(var%mapid)%map)
390  it => maptab(var%mapid)%map(i)
391  j = it%dimno
392  if (j > 0 .and. j <= internal_ndims) then
393  specs(j, 1) = it%start + it%offset
394  specs(j, 2) = it%count
395  if (i > external_ndims) specs(j, 2) = 1
396  specs(j, 3) = it%stride * it%step
397  specs(j, 4) = imap
398  endif
399  imap = imap * it%count
400  enddo
401  end subroutine map_to_internal_specs
402 
403  subroutine map_allocate(map, ndims)
404  ! 次元表エントリに ndims 個のエントリを割り付け初期化する。
405  type(gt_dimmap), pointer:: map(:)
406  integer, intent(in):: ndims
407  if (ndims <= 0) then
408  nullify(map)
409  return
410  endif
411  allocate(map(1:ndims))
412  map(1:ndims)%dimno = -1
413  map(1:ndims)%url = ' '
414  map(1:ndims)%allcount = 0
415  map(1:ndims)%offset = 0
416  map(1:ndims)%step = 1
417  map(1:ndims)%start = 1
418  map(1:ndims)%count = 0
419  map(1:ndims)%stride = 1
420  map(1:ndims)%scalar = .false.
421  end subroutine map_allocate
422 
423  subroutine map_apply(var, map)
424  ! 変数 var にマップ表 map を組み合わせる
425  use gtdata_types, only: gt_variable
426  type(gt_variable), intent(inout):: var
427  type(gt_dimmap), pointer:: map(:)
428  type(gt_dimmap), pointer:: tmpmap(:), varmap
429  integer:: i, nd
430  nd = size(map)
431  allocate(tmpmap(nd))
432  do, i = 1, nd
433  tmpmap(i)%allcount = map(i)%allcount
434  tmpmap(i)%count = map(i)%count
435  if (map(i)%dimno > 0) then
436  varmap => maptab(var%mapid)%map(map(i)%dimno)
437  tmpmap(i)%url = varmap%url
438  tmpmap(i)%dimno = varmap%dimno
439  tmpmap(i)%offset = varmap%offset + map(i)%offset
440  tmpmap(i)%step = varmap%step * map(i)%step
441  else
442  tmpmap(i)%url = map(i)%url
443  tmpmap(i)%dimno = 0
444  tmpmap(i)%offset = map(i)%offset
445  tmpmap(i)%step = map(i)%step
446  endif
447  enddo
448  deallocate(map)
449  map => tmpmap
450  end subroutine map_apply
451 
452  subroutine map_resize(var, ndims)
453  ! 変数 var の次元表の大きさを変える
454  use gtdata_types, only: gt_variable
455  type(gt_variable), intent(in):: var
456  integer, intent(in):: ndims
457  type(gt_dimmap), pointer:: newmap(:)
458  type(gt_dimmap), pointer:: tmpmap(:)
459  integer:: n
460  if (associated(maptab(var%mapid)%map)) then
461  tmpmap => maptab(var%mapid)%map
462  call map_allocate(newmap, ndims)
463  n = min(size(tmpmap), ndims)
464  newmap(1:n) = tmpmap(1:n)
465  deallocate(tmpmap)
466  maptab(var%mapid)%map => newmap
467  newmap(n+1:ndims)%dimno = -1
468  newmap(n+1:ndims)%url = ' '
469  newmap(n+1:ndims)%allcount = 0
470  newmap(n+1:ndims)%offset = 0
471  newmap(n+1:ndims)%step = 1
472  newmap(n+1:ndims)%start = 1
473  newmap(n+1:ndims)%count = 0
474  newmap(n+1:ndims)%stride = 1
475  else
476  call map_allocate(maptab(var%mapid)%map, ndims)
477  n = 1
478  endif
479  end subroutine map_resize
480 
481  subroutine gtvar_dump(var)
482  ! 変数のプロパティを出力
483  use gtdata_types, only: gt_variable
485  use dc_trace, only: debug, dbgmessage
486  type(gt_variable), intent(in):: var
487  integer:: idim, imap
488  logical:: dbg_mode
489  continue
490  call debug( dbg_mode )
491  if (.not. dbg_mode) return
492  imap = var%mapid
493  if (imap < 1 .or. imap > size(maptab)) then
494  call dbgmessage('[gt_variable %d: invalid id]', i=(/imap/))
495  return
496  endif
497  if (associated(maptab(imap)%map)) then
498  call dbgmessage('[gt_variable %d: ndims=%d, map.size=%d]', &
499  & i=(/imap, maptab(imap)%ndims, size(maptab(imap)%map)/))
500  do, idim = 1, size(maptab(imap)%map)
501  call dbgmessage('[dim%d dimno=%d ofs=%d step=%d' &
502  &// ' all=%d start=%d count=%d stride=%d url=%c]', &
503  & c1=trim(maptab(imap)%map(idim)%url), &
504  & i=(/idim, maptab(imap)%map(idim)%dimno, &
505  & maptab(imap)%map(idim)%offset, &
506  & maptab(imap)%map(idim)%step, &
507  & maptab(imap)%map(idim)%allcount, &
508  & maptab(imap)%map(idim)%start, &
509  & maptab(imap)%map(idim)%count, &
510  & maptab(imap)%map(idim)%stride/))
511  enddo
512  else
513  call dbgmessage('[gt_variable %d: ndims=%d, map=null]', &
514  & i=(/imap, maptab(imap)%ndims/))
515  endif
516  call vartable_dump(maptab(imap)%vid)
517  end subroutine gtvar_dump
518 
519  integer function dimord_skip_compact(dimord, map) result(result)
520  ! 次元表の中で非縮退次元だけを数えた次元番号 dimord の次元を
521  ! 特定し、外部向けの次元番号を返す。
522  use dc_trace, only: dbgmessage
523  integer, intent(in):: dimord
524  type(gt_dimmap), intent(in):: map(:)
525  integer:: nd, id
526  result = -1
527  nd = 0
528  do, id = 1, size(map)
529  if (map(id)%count < 2) cycle
530  nd = nd + 1
531  if (nd < dimord) cycle
532  result = id
533  call dbgmessage('compact dim skip: %d <= %d', i=(/result, dimord/))
534  exit
535  enddo
536  end function dimord_skip_compact
537 
538 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
integer, parameter, private maptab_init_size
subroutine map_dup(var, source_var)
type(map_table_entry), dimension(:), allocatable, target, save, private maptab
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)
subroutine, public dbgmessage(fmt, i, r, d, L, n, c1, c2, c3, ca)
Definition: dc_trace.f90:509
integer, parameter, public vtb_class_unused
subroutine, public vartablelookup(vid, class, cid)
種別型パラメタを提供します。
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
integer, parameter, public vtb_class_memory
subroutine map_allocate(map, ndims)
subroutine, public var_class(var, class, cid)
subroutine map_resize(var, ndims)
integer, parameter, public string
文字列を保持する 文字型変数の種別型パラメタ
Definition: dc_types.f90:118