gdncvaropenbydimord.f90
Go to the documentation of this file.
1 !== Open GD_NC_VARIABLE of dimension by dimord
2 !
3 ! Authors:: Yasuhiro MORIKAWA, Eizi TOYODA
4 ! Version:: $Id: gdncvaropenbydimord.f90,v 1.2 2009-05-25 09:51:59 morikawa Exp $
5 ! Tag Name:: $Name: $
6 ! Copyright:: Copyright (C) GFD Dennou Club, 2000-2006. All rights reserved.
7 ! License:: See COPYRIGHT[link:../../COPYRIGHT]
8 !
9 
10 subroutine gdncvaropenbydimord(var, src_var, dimord, err)
11  !
12  ! 既に開かれた変数 *src_var* の *dimord* 番目の次元にあたる変数を
13  ! 開き *var* に格納します。
14  !
15  ! 変数を開く際にエラーが生じた場合、メッセージを出力してプログラムは
16  ! 強制終了します。*err* を与えてある場合にはこの引数に .true.
17  ! が返り、プログラムは終了しません。
18  !
22  use netcdf, only: nf90_max_name, nf90_noerr, nf90_einval, nf90_enotvar, &
23  & nf90_inquire_dimension, nf90_inq_varid
24  use dc_error
25  implicit none
26  type(gd_nc_variable), intent(out):: var
27  type(gd_nc_variable), intent(in):: src_var
28  integer, intent(in):: dimord
29  logical, intent(out), optional:: err
30  type(gd_nc_variable_entry):: src_ent
31  type(gd_nc_variable_search):: ent
32  character(len = NF90_MAX_NAME):: dimname
33  integer:: stat
34 continue
35  stat = vtable_lookup(src_var, src_ent)
36  if (stat /= nf90_noerr) goto 999
37  !
38  if (dimord <= 0) then
39  !
40  ! dimord == 0 の特別扱いは廃止 (gtdata 層でやる)。
41  ! (同じハンドルを返すと参照カウントが増えないので多重クローズが
42  ! 禁止されてしまい、Open のセマンティクスに不適合なため。)
43  !
44  var = src_var
45  stat = nf90_einval
46  goto 999
47  endif
48  if (.not. associated(src_ent%dimids)) then
49  stat = gt_enomoredims
50  goto 999
51  else if (dimord > size(src_ent%dimids)) then
52  stat = gt_enomoredims
53  goto 999
54  endif
55  !
56  ! 決定された次元 id を使って変数 id の確定を試みる
57  !
58  ent%fileid = src_ent%fileid
59  ent%dimid = src_ent%dimids(dimord)
60  stat = nf90_inquire_dimension(ent%fileid, ent%dimid, name = dimname)
61  if (stat /= nf90_noerr) goto 999
62  stat = nf90_inq_varid(ent%fileid, dimname, ent%varid)
63  if (stat == nf90_enotvar) then
64  ! 変数が対応しなければ無視するだけ
65  ent%varid = 0
66  stat = nf90_noerr
67  else if (stat /= nf90_noerr) then
68  goto 999
69  endif
70  !
71  ! 成功しそうなので参照カウント値を増やす
72  call gdncfilereopen(ent%fileid)
73  ! 不適合な変数が開かれている場合は vtable_add が変数部を無視する。
74  stat = vtable_add(var, ent)
75  if (stat /= nf90_noerr) then
76  call gdncfileclose(ent%fileid)
77  goto 999
78  endif
79 
80 999 continue
81  call storeerror(stat, 'GDNcVarOpenByDimOrd', err, cause_i=dimord)
82 end subroutine
integer function, public vtable_lookup(var, entry)
subroutine gdncfileclose(fileid, err)
subroutine gdncfilereopen(fileid, err)
subroutine, public storeerror(number, where, err, cause_c, cause_i)
Definition: dc_error.f90:830
subroutine gdncvaropenbydimord(var, src_var, dimord, err)
integer function, public vtable_add(var, entry)
integer, parameter, public gt_enomoredims
Definition: dc_error.f90:528