subroutine varinfo_init
!==== Dependency
use type_mod, only: STRING, TOKEN, INTKIND, REKIND, DBKIND, NMLARRAY
use nmlfile_mod,only: nmlfile_init, nmlfile_open, nmlfile_close
use dc_types, only: GT_TOKEN => TOKEN, GT_STRING => STRING
use dc_trace, only: DbgMessage, BeginSub, EndSub
use dc_message,only: MessageNotify
use dc_string, only: StrHead, JoinChar
!=end
implicit none
!=begin
!
!==== NAMELIST varinfo_nml
!
!変数に関する基本情報を設定するための NAMELIST で、
!複数の varinfo_nml を用意する事で、複数の変数の設定が可能である。
!変数 varkey はモデル毎に
!内部で設定される物理量のマーカで、具体的には
!((< io_gt4_out_mod >)) モジュールの
!((< io_gt4_out_SetVars >)) サブルーチンでマーカを設定し、
!((< io_gt4_out_Put >)) サブルーチンでデータを file に
!指定されたファイルに出力する。もしも file を与えない、
!または空文字を与えた場合には、((< io_gt4_out_nml >))
!で与えた default_output 変数で指定されたファイルに出力する。
!
!varname, dimnum, dimnames, longname, units, xtype
!は出力される変数に付加される情報である。
!
!StepInterval, OutputStep を与えない、またはゼロ以下の値を
!与えた場合には ((< time_mod >)) の time_nml で与えた
!StepInterval, OutputStep が用いられる。
!
character(STRING) :: varkey = '' ! 変数キー
character(STRING) :: file = '' ! 出力するファイル
character(GT_TOKEN) :: varname = '' ! 変数名
integer(INTKIND) :: dimnum = 0 ! 依存する次元
character(GT_TOKEN) :: dimnames(NMLARRAY) = '' ! 依存する次元
character(GT_STRING) :: longname = '' ! 変数の記述的名称
character(GT_STRING) :: units = '' ! 変数の単位
character(GT_TOKEN) :: xtype = '' ! 変数の型
integer(INTKIND) :: StepInterval = 0 ! 出力ステップ間隔
integer(INTKIND) :: OutputStep = 0 ! 出力回数
namelist /varinfo_nml/ varkey , file , varname , dimnum , dimnames , longname , units , xtype , StepInterval , OutputStep ! 出力ステップ間隔
!
!==== NAMELIST varinfo_attr_nml
!
!変数 varattr の属性情報を与える。
!NAMELIST に複数の varinfo_attr_nml を用意しておく事で
!複数の変数に対し、複数の情報を与える事が可能である。
!与えない場合には属性情報は付加されない。
!
!attrtype には与える属性値の種類を設定する。
!((<URL:http://www.gfd-dennou.org/arch/gtool4/gt4f90io-current/doc/gt_history.htm#derived_gthistoryattr>))
!を参照せよ。なお、arraysize に 1 以上の値を設定すると、
!配列データが優先されて属性値に設定される。
!
character(GT_STRING) :: varattr = '' ! 属性を付加する変数名
character(GT_TOKEN) :: attrname = '' ! 属性名
character(GT_TOKEN) :: attrtype = '' ! 属性値の型
character(GT_STRING) :: cvalue = '' ! 属性の値 (文字)
integer(INTKIND) :: ivalue = 0 ! 属性の値 (整数)
real(REKIND) :: rvalue = 0.0 ! 属性の値 (単精度実数)
real(DBKIND) :: dvalue = 0.0d0 ! 属性の値 (倍精度実数)
logical :: lvalue = .false.! 属性の値 (論理)
integer(INTKIND) :: arraysize= 0 ! 配列のサイズ
integer(INTKIND) :: iarray(NMLARRAY) = 0 ! 属性の値 (整数)
real(REKIND) :: rarray(NMLARRAY) = 0.0 ! 属性の値 (単精度実数)
real(DBKIND) :: darray(NMLARRAY) = 0.0d0! 属性の値 (倍精度実数)
namelist /varinfo_attr_nml/ varattr , attrname , attrtype , cvalue , ivalue , rvalue , dvalue , lvalue , arraysize , iarray , rarray , darray ! 属性の値 (倍精度実数)
!=end
!-------------------------------------------------------------------
! 変数情報の一時格納変数
!-------------------------------------------------------------------
type(VAR_INFO) , allocatable :: vars_store_tmp(:)
! type(GT_HISTORY_ATTR), allocatable :: attrs_tmp(:)
!----------------------------------------------------------------
! 汎用変数
!----------------------------------------------------------------
integer(INTKIND) :: i, j, k
logical :: err
integer(INTKIND) :: nmlstat, nmlunit
logical :: nmlreadable, next
character(TOKEN) :: position
character(STRING) :: var_name
character(STRING), parameter:: subname = "varinfo_init"
continue
!----------------------------------------------------------------
! Check Initialization
!----------------------------------------------------------------
call BeginSub(subname)
if (varinfo_initialized) then
call EndSub( subname, '%c is already called', c1=trim(subname) )
return
else
varinfo_initialized = .true.
endif
!----------------------------------------------------------------
! Version identifier
!----------------------------------------------------------------
call DbgMessage('%c :: %c', c1=trim(version), c2=trim(tagname))
!----------------------------------------------------------------
! read varinfo_nml
!----------------------------------------------------------------
if ( allocated(vars_store) ) then
deallocate(vars_store)
endif
call nmlfile_init
call nmlfile_open(nmlunit, nmlreadable)
if (.not. nmlreadable) then
call DbgMessage('Not Read NAMELIST varinfo_nml')
call MessageNotify('W', subname, 'Can not Read NAMELIST varinfo_nml.')
else
i = 0
j = 0
next = .false.
varinfo_nml_input : do
i = i + 1
call DbgMessage('NAMELIST varinfo_nml Input, <%d> time', i=(/i/))
! 初期化
varkey = '' ! 変数キー
file = '' ! 出力するファイル
varname = '' ! 変数名
dimnum = 0 ! 依存する次元
dimnames(:)= '' ! 依存する次元
longname = '' ! 変数の記述的名称
units = '' ! 変数の単位
xtype = '' ! 変数の型
StepInterval = 0 ! 出力ステップ間隔
OutputStep = 0 ! 出力回数
! read nml
read(nmlunit, nml=varinfo_nml, iostat=nmlstat)
call DbgMessage('Stat of NAMELIST varinfo_nml Input is <%d>', i=(/nmlstat/))
write(0, nml=varinfo_nml)
! Inquire access position
inquire(nmlunit, position=position)
if ( trim(position) /= 'APPEND' .and. nmlstat == 0 ) then
next = .true.
else
next = .false.
endif
! 有効でない値を含むものに関しては無視。
if (varkey == '' ) then
call DbgMessage('var is blank. so this varinfo_nml is ignored.')
if (next) cycle varinfo_nml_input
if (.not. next) exit varinfo_nml_input
elseif (dimnum < 1) then
call DbgMessage('dimnum < 1. so this varinfo_nml is ignored.')
if (next) cycle varinfo_nml_input
if (.not. next) exit varinfo_nml_input
endif
!--------------------------------------------------------------
! vars_store への格納
!--------------------------------------------------------------
! vars_store の準備・拡張
j = j + 1
if ( .not. allocated(vars_store) ) then
allocate( vars_store(1) )
! 初期化
!!$ if ( associated(vars_store(1)%attrs) ) then
!!$ deallocate( vars_store(1)%attrs )
!!$ endif
else
allocate( vars_store_tmp(j-1) )
call varinfo_copy( vars_store(1:j-1), vars_store_tmp(1:j-1) )
deallocate( vars_store )
allocate( vars_store(j) )
call varinfo_copy( vars_store_tmp(1:j-1), vars_store(1:j-1) )
deallocate( vars_store_tmp )
! 初期化
!!$ if ( associated(vars_store(j)%attrs) ) then
!!$ deallocate( vars_store(j)%attrs )
!!$ endif
endif
vars_store(j)%varkey = varkey
vars_store(j)%file = file
vars_store(j)%StepInterval = StepInterval
vars_store(j)%OutputStep = OutputStep
call HistoryVarinfoCreate(vars_store(j)%varinfo, varname, dimnames(1:dimnum), longname, units, xtype)
!!$ vars_store(j)%varinfo%name = varname
!!$ allocate( vars_store(j)%varinfo%dims( dimnum ) )
!!$ vars_store(j)%varinfo%dims = dimnames(1:dimnum)
!!$ vars_store(j)%varinfo%longname = longname
!!$ vars_store(j)%varinfo%units = units
!!$ vars_store(j)%varinfo%xtype = xtype
if (.not. next) exit varinfo_nml_input
next = .false. ! 次回のための初期化
enddo varinfo_nml_input
endif
call nmlfile_close
!----------------------------------------------------------------
! read varinfo_attr_nml
!----------------------------------------------------------------
call nmlfile_init
call nmlfile_open(nmlunit, nmlreadable)
if (.not. nmlreadable) then
call DbgMessage('Not Read NAMELIST varinfo_attr_nml')
call MessageNotify('W', subname, 'Can not Read NAMELIST varinfo_attr_nml.')
else
i = 0
next = .false.
varinfo_attr_nml_input : do
if ( .not. allocated(vars_store) ) then
call DbgMessage('variables are not defined, so varinfo_attr_nml is ignored')
exit varinfo_attr_nml_input
endif
i = i + 1
call DbgMessage('NAMELIST varinfo_attr_nml Input, <%d> time', i=(/i/))
! 初期化
varattr = '' ! 属性を付加する変数名
attrname = '' ! 属性名
attrtype = '' ! 属性値の型
cvalue = '' ! 属性の値 (文字)
ivalue = 0 ! 属性の値 (整数)
rvalue = 0.0 ! 属性の値 (単精度実数)
dvalue = 0.0d0 ! 属性の値 (倍精度実数)
lvalue = .false.! 属性の値 (論理)
arraysize = 0 ! 配列のサイズ
iarray(:) = 0 ! 属性の値 (整数)
rarray(:) = 0.0 ! 属性の値 (単精度実数)
darray(:) = 0.0d0 ! 属性の値 (倍精度実数)
! read nml
read(nmlunit, nml=varinfo_attr_nml, iostat=nmlstat)
call DbgMessage('Stat of NAMELIST varinfo_attr_nml Input is <%d>', i=(/nmlstat/))
write(0, nml=varinfo_attr_nml)
! Inquire access position
inquire(nmlunit, position=position)
if ( trim(position) /= 'APPEND' .and. nmlstat == 0 ) then
next = .true.
else
next = .false.
endif
! 有効でない値を含むものに関しては無視。
if (varattr == '') then
call DbgMessage('varattr is blank. so this varinfo_attr_nml is ignored.')
if (next) cycle
if (.not. next) exit varinfo_attr_nml_input
elseif (attrname == '') then
call DbgMessage('attrname is blank. so this varinfo_attr_nml is ignored.')
if (next) cycle
if (.not. next) exit varinfo_attr_nml_input
elseif (attrtype == '') then
call DbgMessage('attrtype is blank. so this varinfo_attr_nml is ignored.')
if (next) cycle
if (.not. next) exit varinfo_attr_nml_input
endif
! varattr を vars_store のどこに格納すべきか探査して j に格納。
! 対応するものが vars_store 内に無い場合は無視。
j = 0
do j = 1, size(vars_store)
err = .false.
call HistoryVarinfoInquire(vars_store(j)%varinfo, name=var_name)
call DbgMessage('Search varname=<%c> in vars_store (<%c>)', c1=trim(varattr) , c2=trim(var_name) )
if ( trim(varattr) == trim(var_name) ) then
exit
endif
err = .true.
enddo
if (err) then
call DbgMessage('variable <%c> is not defined in varinfo_nml.'// ' So this varinfo_attr_nml is ignored.' , c1=trim(varattr) )
if (next) cycle
if (.not. next) exit
endif
!--------------------------------------------------------------
! vars_store%attrs への格納
!--------------------------------------------------------------
! attrs(:) の拡張
!!$ if ( .not. associated(vars_store(j)%attrs) ) then
!!$ allocate( vars_store(j)%attrs(1) )
!!$ k = 1
!!$ else
!!$ k = size( vars_store(j)%attrs ) + 1
!!$ allocate( attrs_tmp(k-1) )
!!$ call varinfo_attrs_copy(vars_store(j)%attrs(1:k-1), attrs_tmp(1:k-1))
!!$ deallocate( vars_store(j)%attrs )
!!$ allocate( vars_store(j)%attrs(k) )
!!$ call varinfo_attrs_copy(attrs_tmp(1:k-1), vars_store(j)%attrs(1:k-1))
!!$ deallocate( attrs_tmp )
!!$ endif
!!$
!!$ if (arraysize > 0) then
!!$ call varinfo_attrs_init(vars_store(j)%attrs(k))
!!$
!!$ deallocate( vars_store(j)%attrs(k)%iarray )
!!$ deallocate( vars_store(j)%attrs(k)%rarray )
!!$ deallocate( vars_store(j)%attrs(k)%darray )
!!$
!!$ allocate( vars_store(j)%attrs(k)%iarray( arraysize ) )
!!$ allocate( vars_store(j)%attrs(k)%rarray( arraysize ) )
!!$ allocate( vars_store(j)%attrs(k)%darray( arraysize ) )
!!$
!!$ vars_store(j)%attrs(k)%array = .true.
!!$
!!$ else
!!$ call varinfo_attrs_init(vars_store(j)%attrs(k))
!!$ endif
if (StrHead(attrtype, 'char')) then
call HistoryVarinfoAddAttr(vars_store(j) % varinfo, attrname, cvalue)
elseif (StrHead(attrtype, 'logical')) then
call HistoryVarinfoAddAttr(vars_store(j) % varinfo, attrname, lvalue)
elseif (StrHead(attrtype, 'int')) then
if (arraysize > 0) then
call HistoryVarinfoAddAttr(vars_store(j) % varinfo, attrname, iarray)
else
call HistoryVarinfoAddAttr(vars_store(j) % varinfo, attrname, ivalue)
end if
elseif (StrHead(attrtype, 'real')) then
if (arraysize > 0) then
call HistoryVarinfoAddAttr(vars_store(j) % varinfo, attrname, rarray)
else
call HistoryVarinfoAddAttr(vars_store(j) % varinfo, attrname, rvalue)
end if
elseif (StrHead(attrtype, 'double')) then
if (arraysize > 0) then
call HistoryVarinfoAddAttr(vars_store(j) % varinfo, attrname, darray)
else
call HistoryVarinfoAddAttr(vars_store(j) % varinfo, attrname, dvalue)
end if
end if
!!$ vars_store(j)%attrs(k)%attrname = attrname
!!$ vars_store(j)%attrs(k)%attrtype = attrtype
!!$ vars_store(j)%attrs(k)%cvalue = cvalue
!!$ vars_store(j)%attrs(k)%ivalue = ivalue
!!$ vars_store(j)%attrs(k)%rvalue = rvalue
!!$ vars_store(j)%attrs(k)%dvalue = dvalue
!!$ vars_store(j)%attrs(k)%lvalue = lvalue
!!$ vars_store(j)%attrs(k)%iarray(1:max(1,arraysize)) = iarray(1:max(1,arraysize))
!!$ vars_store(j)%attrs(k)%rarray(1:max(1,arraysize)) = rarray(1:max(1,arraysize))
!!$ vars_store(j)%attrs(k)%darray(1:max(1,arraysize)) = darray(1:max(1,arraysize))
!!$
if (.not. next) exit varinfo_attr_nml_input
next = .false. ! 次回のための初期化
enddo varinfo_attr_nml_input
end if
call nmlfile_close
call EndSub( subname )
end subroutine varinfo_init