| Class | axis_y_mod |
| In: |
shared/axis_y.f90
|
| Subroutine : |
subroutine axis_y_end()
!==== Dependency
use dc_trace, only: DbgMessage, BeginSub, EndSub
use gt4_history, only: HistoryAxisClear
!=end
implicit none
character(STRING), parameter:: subname = "axis_y_end"
continue
!----------------------------------------------------------------
! Check Initialization
!----------------------------------------------------------------
call BeginSub(subname)
if ( .not. axis_y_initialized) then
call EndSub( subname, 'axis_y_init was not called', c1=trim(subname) )
return
else
axis_y_initialized = .false.
axis_y_data_from_spectral = .false.
axis_y_data_from_manual = .false.
axis_y_data_from_netcdf = ''
!--------------------------------------------------------------
! Initialize y_Dim
!--------------------------------------------------------------
y_Dim%stored = .false.
call HistoryAxisClear(y_Dim % axisinfo)
!!$ y_Dim%axisinfo%name = ''
!!$ y_Dim%axisinfo%length = 0
!!$ y_Dim%axisinfo%longname = ''
!!$ y_Dim%axisinfo%units = ''
!!$ y_Dim%axisinfo%xtype = ''
!!$ if ( associated(y_Dim%a_Dim) ) then
!!$ deallocate( y_Dim%a_Dim )
!!$ endif
!!$
!!$ if ( associated(y_Dim%attrs) ) then
!!$ deallocate( y_Dim%attrs )
!!$ end if
!--------------------------------------------------------------
! Initialize y_Dim_Weight
!--------------------------------------------------------------
call HistoryAxisClear(y_Dim_Weight % axisinfo)
!!$ y_Dim_Weight%axisinfo%name = ''
!!$ y_Dim_Weight%axisinfo%length = 0
!!$ y_Dim_Weight%axisinfo%longname = ''
!!$ y_Dim_Weight%axisinfo%units = ''
!!$ y_Dim_Weight%axisinfo%xtype = ''
!!$
!!$ if ( associated(y_Dim%a_Dim) ) then
!!$ deallocate( y_Dim%a_Dim )
!!$ end if
!!$
!!$ if ( associated(y_Dim%attrs) ) then
!!$ deallocate( y_Dim%attrs )
!!$ end if
endif
call EndSub( subname )
end subroutine axis_y_end
| Subroutine : |
This procedure input/output NAMELIST#axis_y_nml, NAMELIST#axis_y_attr_nml .
subroutine axis_y_init
!==== Dependency
use type_mod, only: STRING, TOKEN, INTKIND, REKIND, DBKIND, NMLARRAY
use nmlfile_mod,only: nmlfile_init, nmlfile_open, nmlfile_close
use grid_3d_mod,only: jm, grid_3d_init
use spml_mod, only: spml_init
! use axis_type_mod, only : axis_attrs_copy, axis_attrs_init
! use gt4_history, only: GT_HISTORY_ATTR
use dc_types, only: GT_TOKEN => TOKEN, GT_STRING => STRING
use dc_url, only: GT_ATMARK, GT_QUESTION
use dc_trace, only: DbgMessage, BeginSub, EndSub
use dc_message,only: MessageNotify
use gt4_history, only: HistoryAxisCreate, HistoryAxisAddAttr
use dc_string, only: StrHead
!=end
implicit none
!-------------------------------------------------------------------
! 変数定義
!-------------------------------------------------------------------
!=begin
!
!==== NAMELIST
!
!X 軸の次元変数に関する情報を与える。
!値を与えないものに関しては以下のデフォルトの値が用いられる。
!
!変数 decision には Y 軸のデータをどのように与えるかを指定する。
!
!* (({ 'manual' }))
! * Data 配列に格納したデータをそのまま Y 軸として与える。
!
!* (({ 'spectral' }))
! * スペクトル法を用いる事を想定し、 length に応じて自動的に
! Y 軸のデータが決まる。
!
!* gtool4 変数 (例えば (({ 'foo.nc@lon' })) など)
! * 該当する変数から Y 軸のデータを取得する。
!
!変数 length には、 ((<grid_3d_mod>)) の公開要素 ((< jm >)) と同じ
!値を与えなければならない。
!
character(TOKEN) :: name = 'lon' ! 次元変数名
integer(INTKIND) :: length = 64 ! 次元長 (配列サイズ)
character(STRING) :: longname = 'Longitude' ! 次元変数の記述的名称
character(STRING) :: units = 'degrees_east' ! 次元変数の単位
character(TOKEN) :: xtype = 'float' ! 次元変数の型
character(STRING) :: decision = 'spectral' ! 次元データの取得方法
real(REKIND) :: Data(NMLARRAY) = 0.0 ! 次元データ入力用
namelist /axis_y_nml/ name , length , longname , units , xtype , decision , Data ! 次元データ
!=end
!=begin
!
!X 軸の次元変数の属性に関する情報を与える。
!NAMELIST に複数の axis_y_attr_nml を用意しておく事で
!複数の情報を与える事が可能である。
!与えない場合には属性情報は付加されない。
!
!attrtype には与える属性値の種類を設定する。
!((<URL:http://www.gfd-dennou.org/arch/gtool4/gt4f90io-current/doc/gt_history.htm#derived_gthistoryattr>))
!を参照せよ。なお、arraysize に 1 以上の値を設定すると、
!配列データが優先されて属性値に設定される。
!
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 /axis_y_attr_nml/ attrname , attrtype , cvalue , ivalue , rvalue , dvalue , lvalue , arraysize , iarray , rarray , darray ! 属性の値 (倍精度実数)
!=end
!-----------------------------------------------------------------
! 変数情報の一時格納変数
!-----------------------------------------------------------------
! type(GT_HISTORY_ATTR), allocatable :: attrs_tmp(:)
!-----------------------------------------------------------------
! 汎用変数
!-----------------------------------------------------------------
integer(INTKIND) :: i, k
integer(INTKIND) :: nmlstat, nmlunit
logical :: nmlreadable, next
character(TOKEN) :: position
character(STRING), parameter:: subname = "axis_y_init"
continue
!----------------------------------------------------------------
! Check Initialization
!----------------------------------------------------------------
call BeginSub(subname)
if (axis_y_initialized) then
call EndSub( subname, '%c is already called', c1=trim(subname) )
return
else
axis_y_initialized = .true.
endif
!----------------------------------------------------------------
! Version identifier
!----------------------------------------------------------------
call DbgMessage('%c :: %c', c1=trim(version), c2=trim(tagname))
!-----------------------------------------------------------------
! Initialize Dependent Modules
!-----------------------------------------------------------------
call grid_3d_init
call spml_init
!----------------------------------------------------------------
! read axis_y_nml
!----------------------------------------------------------------
call nmlfile_init
call nmlfile_open(nmlunit, nmlreadable)
if (nmlreadable) then
read(nmlunit, nml=axis_y_nml, iostat=nmlstat)
call DbgMessage('Stat of NAMELIST axis_y_nml Input is <%d>', i=(/nmlstat/))
write(0, nml=axis_y_nml)
else
call DbgMessage('Not Read NAMELIST axis_y_nml')
call MessageNotify('W', subname, 'Can not Read NAMELIST axis_y_nml. Force Use Default Value.')
end if
call nmlfile_close
y_Dim%stored = .false.
! 次元変数の情報を構造型変数 y_Dim への代入
call HistoryAxisCreate(y_Dim % axisinfo, name, length, longname, units, xtype)
!!$ y_Dim%axisinfo%name = name
!!$ y_Dim%axisinfo%length = length
!!$ y_Dim%axisinfo%longname = longname
!!$ y_Dim%axisinfo%units = units
!!$ y_Dim%axisinfo%xtype = xtype
!!$ allocate( y_Dim%a_Dim(y_Dim%axisinfo%length) )
allocate( y_Dim%a_Dim(length) )
! 次元変数の情報を構造型変数 y_Dim への代入
select case(decision)
! manual: NAMELIST の Data で入力
case('manual')
y_Dim%a_Dim(:) = 0
y_Dim%a_Dim(1:length) = Data(1:length)
y_Dim%stored = .true.
axis_y_data_from_manual = .true.
! spectral: SPMODEL の wa_module により生成
case('spectral')
axis_y_data_from_spectral = .true.
y_Dim%stored = .false.
! foo.nc@lon: foo.nc ファイルの lon 変数から取得
! その他 : spectral と同じに
case default
! 文字の中に '@' か '?' が含まれる場合は gtool4 変数として
! 認識し、その nc ファイルから変数情報をコピーする。
if ( index(decision, GT_ATMARK) > 0 .or. index(decision, GT_QUESTION) > 0) then
axis_y_data_from_netcdf = decision
y_Dim%stored = .false.
! それ以外は 'spectral' と同じように処理
else
axis_y_data_from_spectral = .true.
y_Dim%stored = .false.
endif
end select
!----------------------------------------------------------------
! read axis_y_attr_nml
!----------------------------------------------------------------
call nmlfile_init
call nmlfile_open(nmlunit, nmlreadable)
if (.not. nmlreadable) then
call DbgMessage('Not Read NAMELIST axis_y_attr_nml')
call MessageNotify('W', subname, 'Can not Read NAMELIST axis_y_attr_nml.')
else
i = 0
next = .false.
axis_y_attr_nml_input : do
i = i + 1
call DbgMessage('NAMELIST axis_y_attr_nml Input, <%d> time', i=(/i/))
! 初期化
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=axis_y_attr_nml, iostat=nmlstat)
call DbgMessage('Stat of NAMELIST axis_y_attr_nml Input is <%d>', i=(/nmlstat/))
write(0, nml=axis_y_attr_nml)
! Inquire access position
inquire(nmlunit, position=position)
if ( trim(position) /= 'APPEND' .and. nmlstat == 0 ) then
next = .true.
else
next = .false.
endif
! 有効でない値を含むものに関しては無視。
if (attrname == '') then
call DbgMessage('attrname is blank. so this axis_y_attr_nml is ignored.')
if (next) cycle
if (.not. next) exit
elseif (attrtype == '') then
call DbgMessage('attrtype is blank. so this axis_y_attr_nml is ignored.')
if (next) cycle
if (.not. next) exit
endif
!-----------------------------------------------------------
! y_Dim%attrs への格納
!-----------------------------------------------------------
! attrs(:) の拡張
!!$ if ( .not. associated(y_Dim%attrs) ) then
!!$ allocate( y_Dim%attrs(1) )
!!$ k = 1
!!$ else
!!$ k = size( y_Dim%attrs ) + 1
!!$ ! 配列データの領域確保
!!$ allocate( attrs_tmp(k-1) )
!!$ call axis_attrs_copy(from=y_Dim%attrs(1:k-1), to=attrs_tmp(1:k-1))
!!$ deallocate( y_Dim%attrs )
!!$ allocate( y_Dim%attrs(k) )
!!$ call axis_attrs_copy(from=attrs_tmp(1:k-1), to=y_Dim%attrs(1:k-1))
!!$ deallocate( attrs_tmp )
!!$ endif
!!$
!!$ if (arraysize > 0) then
!!$ call axis_attrs_init(y_Dim%attrs(k))
!!$
!!$ deallocate( y_Dim%attrs(k)%iarray )
!!$ deallocate( y_Dim%attrs(k)%rarray )
!!$ deallocate( y_Dim%attrs(k)%darray )
!!$
!!$ allocate( y_Dim%attrs(k)%iarray( arraysize ) )
!!$ allocate( y_Dim%attrs(k)%rarray( arraysize ) )
!!$ allocate( y_Dim%attrs(k)%darray( arraysize ) )
!!$
!!$ y_Dim%attrs(k)%array = .true.
!!$
!!$ else
!!$ call axis_attrs_init(y_Dim%attrs(k))
!!$ endif
if (StrHead(attrtype, 'char')) then
call HistoryAxisAddAttr(y_Dim % axisinfo, attrname, cvalue)
elseif (StrHead(attrtype, 'logical')) then
call HistoryAxisAddAttr(y_Dim % axisinfo, attrname, lvalue)
elseif (StrHead(attrtype, 'int')) then
if (arraysize > 0) then
call HistoryAxisAddAttr(y_Dim % axisinfo, attrname, iarray)
else
call HistoryAxisAddAttr(y_Dim % axisinfo, attrname, ivalue)
end if
elseif (StrHead(attrtype, 'real')) then
if (arraysize > 0) then
call HistoryAxisAddAttr(y_Dim % axisinfo, attrname, rarray)
else
call HistoryAxisAddAttr(y_Dim % axisinfo, attrname, rvalue)
end if
elseif (StrHead(attrtype, 'double')) then
if (arraysize > 0) then
call HistoryAxisAddAttr(y_Dim % axisinfo, attrname, darray)
else
call HistoryAxisAddAttr(y_Dim % axisinfo, attrname, dvalue)
end if
end if
!!$ y_Dim%attrs(k)%attrname = attrname
!!$ y_Dim%attrs(k)%attrtype = attrtype
!!$ y_Dim%attrs(k)%cvalue = cvalue
!!$ y_Dim%attrs(k)%ivalue = ivalue
!!$ y_Dim%attrs(k)%rvalue = rvalue
!!$ y_Dim%attrs(k)%dvalue = dvalue
!!$ y_Dim%attrs(k)%lvalue = lvalue
!!$ y_Dim%attrs(k)%iarray(1:max(1,arraysize)) = iarray(1:max(1,arraysize))
!!$ y_Dim%attrs(k)%rarray(1:max(1,arraysize)) = rarray(1:max(1,arraysize))
!!$ y_Dim%attrs(k)%darray(1:max(1,arraysize)) = darray(1:max(1,arraysize))
!!$ call DbgMessage('y_Dim-attrs(%d) [attrname=<%c> ' // &
!!$ & 'attrtype=<%c> array=<%b> cvalue=<%c> ' // &
!!$ & 'ivalue=<%d> rvalue=<%r> dvalue=<%f> ' // &
!!$ & 'iarray(1:%d)=<%d, ...> ' // &
!!$ & 'rarray(1:%d)=<%r, ...> darray(1:%d)=<%f, ...>' , &
!!$ & c1=trim( y_Dim%attrs(k)%attrname ) , &
!!$ & c2=trim( y_Dim%attrs(k)%attrtype ) , &
!!$ & c3=trim( y_Dim%attrs(k)%cvalue ) , &
!!$ & i=(/ k, y_Dim%attrs(k)%ivalue , &
!!$ & size(y_Dim%attrs(k)%iarray) , &
!!$ & y_Dim%attrs(k)%iarray , &
!!$ & size(y_Dim%attrs(k)%rarray) , &
!!$ & size(y_Dim%attrs(k)%darray) &
!!$ & /) , &
!!$ & r=(/y_Dim%attrs(k)%rvalue, y_Dim%attrs(k)%rarray/) , &
!!$ & d=(/y_Dim%attrs(k)%dvalue, y_Dim%attrs(k)%darray/) , &
!!$ & l=(/y_Dim%attrs(k)%lvalue/) )
if (.not. next) exit axis_y_attr_nml_input
next = .false. ! 次回のための初期化
enddo axis_y_attr_nml_input
end if
call nmlfile_close
!----------------------------------------------------------------
! grid_3d_mod との整合性をチェック
!----------------------------------------------------------------
call grid_3d_init
if (length /= jm) then
call MessageNotify('E', subname, message='axis length is inconsistent with im in grid_3d_mod')
endif
!----------------------------------------------------------------
! 例外処理
!----------------------------------------------------------------
if (length < 1) then
call MessageNotify('E', subname, message='Invalid grid number.')
endif
call EndSub( subname )
end subroutine axis_y_init
| Subroutine : | |||
| Dim : | type(AXISINFO), intent(inout)
|
subroutine axis_y_manual(Dim)
!==== Dependency
use axis_type_mod, only: axis_type_copy
use spml_mod, only: wa_module_y_Lat => y_Lat
use dc_trace, only: DbgMessage, BeginSub, EndSub
!=end
implicit none
!=begin
!==== In/Out
!
type(AXISINFO), intent(inout) :: Dim ! 次元情報を包括する変数
!=end
character(STRING), parameter:: subname = "axis_y_manual"
continue
!----------------------------------------------------------------
! Check Initialization
!----------------------------------------------------------------
call BeginSub( subname )
if (.not. axis_y_initialized) then
call EndSub( subname, 'Call axis_y_init before call %c', c1=trim(subname) )
return
endif
!----------------------------------------------------------------
! decision が 'manual' でない場合は停止して返す。
!----------------------------------------------------------------
if (.not. axis_y_data_from_manual) then
call EndSub( subname, 'In axis_y_nml, configurated to not generate from manual')
return
endif
!----------------------------------------------------------------
! intent(out) の引数にデータを渡す。
!----------------------------------------------------------------
call axis_type_copy(y_Dim, Dim)
call EndSub( subname )
end subroutine axis_y_manual
| Subroutine : | |||
| Dim : | type(AXISINFO), intent(inout)
|
subroutine axis_y_netcdf(Dim)
!==== Dependency
use axis_type_mod, only: axis_type_copy
use gt4_history,only: HistoryGet
use dc_url , only: UrlSplit
use dc_trace, only: DbgMessage, BeginSub, EndSub
!=end
implicit none
!=begin
!==== In/Out
!
type(AXISINFO), intent(inout) :: Dim ! 次元情報を包括する変数
!=end
character(STRING) :: file
character(STRING) :: varname
character(STRING), parameter:: subname = "axis_y_netcdf"
continue
!----------------------------------------------------------------
! Check Initialization
!----------------------------------------------------------------
call BeginSub( subname )
if (.not. axis_y_initialized) then
call EndSub( subname, 'Call axis_y_init before call %c', c1=trim(subname) )
return
endif
!----------------------------------------------------------------
! decision が '***@lon' などでない場合は停止して返す。
!----------------------------------------------------------------
if (axis_y_data_from_netcdf == '') then
call EndSub( subname, 'In axis_y_nml, configurated to not generate from netcdf data')
return
endif
!----------------------------------------------------------------
! HistoryGet によるデータ入力
!----------------------------------------------------------------
call UrlSplit(axis_y_data_from_netcdf, file=file, var=varname)
call HistoryGet(file, varname, y_Dim%a_Dim)
y_Dim%stored = .true.
! 入力するデータの units が 'degree' か 'radian' か調べ、
! 一時、'radian' に変換し、その後、y_Dim%axisinfo%units が
! 'degree' か 'radian' かに応じて 180. / pi を掛けるかどうか
! 判定すべき
!----------------------------------------------------------------
! intent(out) の引数にデータを渡す。
!----------------------------------------------------------------
call axis_type_copy(y_Dim, Dim)
call EndSub( subname )
end subroutine axis_y_netcdf
| Subroutine : | |||
| Dim : | type(AXISINFO), intent(inout)
|
subroutine axis_y_spectral(Dim)
!==== Dependency
use axis_type_mod, only: axis_type_copy
use constants_mod, only: constants_init, pi
use spml_mod, only: wa_module_y_Lat => y_Lat
use dc_string, only: toChar, StrHead, LChar
use dc_trace, only: DbgMessage, BeginSub, EndSub
use gt4_history, only: HistoryAxisInquire
!=end
implicit none
!=begin
!==== In/Out
!
type(AXISINFO), intent(inout) :: Dim ! 次元情報を包括する変数
!=end
character(STRING) :: units
character(STRING), parameter:: subname = "axis_y_spectral"
continue
!----------------------------------------------------------------
! Check Initialization
!----------------------------------------------------------------
call BeginSub( subname )
if (.not. axis_y_initialized) then
call EndSub( subname, 'Call axis_y_init before call %c', c1=trim(subname) )
return
endif
!----------------------------------------------------------------
! decision が 'spectral' でない場合は停止して返す。
!----------------------------------------------------------------
if (.not. axis_y_data_from_spectral) then
call EndSub( subname, 'In axis_y_nml, configurated to not generate from spmodel')
return
endif
!----------------------------------------------------------------
! Initialize Dependent module
!----------------------------------------------------------------
call constants_init()
!----------------------------------------------------------------
! wa_module からスペクトルデータを受け取る。
!----------------------------------------------------------------
call HistoryAxisInquire(y_Dim % axisinfo, units=units)
if ( StrHead( 'radians', trim(LChar(units)) ) .or. StrHead( 'rad.', trim(LChar(units)) ) ) then
! radian でそのまま代入
y_Dim%a_Dim(:) = wa_module_y_Lat(:)
else
! radian を degree に変換して代入
y_Dim%a_Dim(:) = wa_module_y_Lat(:) * 180. / pi
end if
y_Dim%stored = .true.
call DbgMessage('y_Lat=<%c>', c1=trim(toChar(y_Dim%a_Dim)) )
!----------------------------------------------------------------
! intent(out) の引数にデータを渡す。
!----------------------------------------------------------------
call axis_type_copy(y_Dim, Dim)
call EndSub( subname )
end subroutine axis_y_spectral
| Subroutine : | |||
| Dim_Weight : | type(AXISINFO), intent(out)
|
subroutine axis_y_weight(Dim_Weight)
!==== Dependency
use constants_mod, only: constants_init, pi
use axis_type_mod, only: axis_type_copy
! use axis_type_mod, only: axis_type_copy, axis_attrs_copy, axis_attrs_init
use spml_mod, only: wa_module_y_Lat_Weight => y_Lat_Weight
use grid_3d_mod,only: im
! use gt4_history,only: GT_HISTORY_ATTR
use dc_trace, only: DbgMessage, BeginSub, EndSub
use gt4_history,only: HistoryAxisCopy, HistoryAxisInquire
!=end
implicit none
!=begin
!==== Output
!
type(AXISINFO), intent(out) :: Dim_Weight ! 次元情報を包括する変数
!=end
! type(GT_HISTORY_ATTR), allocatable :: attrs_tmp(:) ! 次元情報を一時格納
character(STRING) :: name, longname
character(STRING), parameter:: subname = "axis_y_weight"
continue
!----------------------------------------------------------------
! Check Initialization
!----------------------------------------------------------------
call BeginSub( subname )
if (.not. axis_y_initialized) then
call EndSub( subname, 'Call axis_y_init before call %c', c1=trim(subname) )
return
endif
!----------------------------------------------------------------
! decision が 'spectral' でない場合は停止して返す。
!----------------------------------------------------------------
if (.not. axis_y_data_from_spectral) then
call EndSub( subname, 'In axis_y_nml, configurated to not generate from spmodel')
return
endif
!----------------------------------------------------------------
! Initialize Dependent module
!----------------------------------------------------------------
call constants_init()
!----------------------------------------------------------------
! y_Dim の情報から y_Dim_Weight の情報を生成。
!----------------------------------------------------------------
call HistoryAxisInquire(y_Dim % axisinfo, name=name, longname=longname)
call HistoryAxisCopy( axis_dest=y_Dim_Weight % axisinfo, axis_src=y_Dim % axisinfo, name=trim(name) // '_weight', longname=trim(longname) // ' weight')
!!$ y_Dim_Weight%axisinfo%name = trim(y_Dim%axisinfo%name) // '_weight'
!!$ y_Dim_Weight%axisinfo%length = y_Dim%axisinfo%length
!!$ y_Dim_Weight%axisinfo%longname = trim(y_Dim%axisinfo%longname) // ' weight'
!!$ y_Dim_Weight%axisinfo%units = trim(y_Dim%axisinfo%units)
!!$ y_Dim_Weight%axisinfo%xtype = trim(y_Dim%axisinfo%xtype)
!!$ if ( associated(y_Dim_Weight%attrs) ) then
!!$ deallocate( y_Dim_Weight%attrs )
!!$ endif
!!$
!!$ ! y_Dim にも情報付加
!!$ if ( associated(y_Dim%attrs) ) then
!!$ allocate( attrs_tmp(size(y_Dim%attrs)) )
!!$ call axis_attrs_copy(from=y_Dim%attrs(:), to=attrs_tmp(:))
!!$ deallocate( y_Dim%attrs )
!!$ allocate( y_Dim%attrs(size(attrs_tmp)+1) )
!!$ call axis_attrs_copy(from=attrs_tmp(:), to=y_Dim%attrs(1:size(attrs_tmp)) )
!!$
!!$ call axis_attrs_init( y_Dim%attrs(size(attrs_tmp)+1) )
!!$ y_Dim%attrs(size(attrs_tmp)+1)%attrname = 'gt_calc_weight'
!!$ y_Dim%attrs(size(attrs_tmp)+1)%attrtype = 'char'
!!$ y_Dim%attrs(size(attrs_tmp)+1)%cvalue = y_Dim_Weight%axisinfo%name
!!$ y_Dim%attrs(size(attrs_tmp)+1)%array = .false.
!!$ else
!!$ allocate( y_Dim%attrs(1) )
!!$ call axis_attrs_init( y_Dim%attrs(1) )
!!$ y_Dim%attrs(1)%attrname = 'gt_calc_weight'
!!$ y_Dim%attrs(1)%attrtype = 'char'
!!$ y_Dim%attrs(1)%cvalue = y_Dim_Weight%axisinfo%name
!!$ y_Dim%attrs(1)%array = .false.
!!$ endif
!----------------------------------------------------------------
! wa_module からスペクトルデータを受け取る。
!----------------------------------------------------------------
allocate( y_Dim_Weight%a_Dim(im) )
! ラジアンを度数に変換して代入
y_Dim_Weight%a_Dim(:) = wa_module_y_Lat_Weight(:) * 180. / pi
y_Dim_Weight%stored = .true.
!----------------------------------------------------------------
! intent(out) の引数にデータを渡す。
!----------------------------------------------------------------
call axis_type_copy(y_Dim_Weight, Dim_Weight)
call EndSub( subname )
end subroutine axis_y_weight