Class | axis_x_mod |
In: |
axis/axis_x.f90
|
Copyright (C) GFD Dennou Club, 2005. All rights reserved.
!=begin
* Developers: Morikawa Yasuhiro * Version: $Id: axis_x.f90,v 1.8 2005/01/19 08:52:24 morikawa Exp $ * Tag Name: $Name: $ * Change History:
This module set axis X or axis Longitude. X 軸または経度軸を設定する。
netCDF データから X 軸を入力する axis_x_netcdf にて、 元データが radians でも degrees でも、そのまま入力されるように なっている。本来は元データの units から判定すべき。
!=end
((<axis_x_init>)) で設定された値を破棄し、デフォルトに戻す。
subroutine axis_x_end () !==== Dependency use dc_trace, only: DbgMessage, BeginSub, EndSub !=end implicit none character(STRING), parameter:: subname = "axis_x_end" continue !---------------------------------------------------------------- ! Check Initialization !---------------------------------------------------------------- call BeginSub(subname) if ( .not. axis_x_initialized) then call EndSub( subname, 'axis_x_init was not called', & & c1=trim(subname) ) return else axis_x_initialized = .false. axis_x_data_from_spectral = .false. axis_x_data_from_manual = .false. axis_x_data_from_netcdf = '' !-------------------------------------------------------------- ! Initialize x_Dim !-------------------------------------------------------------- x_Dim%stored = .false. x_Dim%axisinfo%name = '' x_Dim%axisinfo%length = 0 x_Dim%axisinfo%longname = '' x_Dim%axisinfo%units = '' x_Dim%axisinfo%xtype = '' if ( allocated(x_Dim%a_Dim) ) then deallocate( x_Dim%a_Dim ) endif if ( allocated(x_Dim%attrs) ) then deallocate( x_Dim%attrs ) end if !-------------------------------------------------------------- ! Initialize x_Dim_Weight !-------------------------------------------------------------- x_Dim_Weight%axisinfo%name = '' x_Dim_Weight%axisinfo%length = 0 x_Dim_Weight%axisinfo%longname = '' x_Dim_Weight%axisinfo%units = '' x_Dim_Weight%axisinfo%xtype = '' if ( allocated(x_Dim%a_Dim) ) then deallocate( x_Dim%a_Dim ) end if if ( allocated(x_Dim%attrs) ) then deallocate( x_Dim%attrs ) end if endif call EndSub( subname ) end subroutine
NAMELIST から代入されたデータを X 軸データとして返す。
((< axis_x_init >)) の NAMELIST axis_x_nml の decision 変数で (({ ‘manual’ })) 以外が与えられた場合は値を代入しないで返す。
subroutine axis_x_manual (Dim) !==== Dependency use axis_type_mod, only: axis_type_copy use spml_mod, only: wa_module_x_Lon => x_Lon use dc_trace, only: DbgMessage, BeginSub, EndSub !=end implicit none !=begin !==== In/Out ! type(AXISINFO), intent(inout) :: Dim ! 次元情報を包括する変数 !=end character(STRING), parameter:: subname = "axis_x_manual" continue !---------------------------------------------------------------- ! Check Initialization !---------------------------------------------------------------- call BeginSub( subname ) if (.not. axis_x_initialized) then call EndSub( subname, 'Call axis_x_init before call %c', & & c1=trim(subname) ) return endif !---------------------------------------------------------------- ! decision が 'manual' でない場合は停止して返す。 !---------------------------------------------------------------- if (.not. axis_x_data_from_manual) then call EndSub( subname, & & 'In axis_x_nml, configurated to not generate from manual') return endif !---------------------------------------------------------------- ! intent(out) の引数にデータを渡す。 !---------------------------------------------------------------- call axis_type_copy(x_Dim, Dim) call EndSub( subname ) end subroutine
netCDF データから取得したデータを X 軸のデータとして返す。
現在、取得先のデータの単位に関わらず、そのままデータが 入力される。((<Known Bugs>)) 参照。
((< axis_x_init >)) の NAMELIST axis_x_nml の decision 変数で gtool4 変数以外が与えられた場合は値を代入しないで返す。
subroutine axis_x_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_x_netcdf" continue !---------------------------------------------------------------- ! Check Initialization !---------------------------------------------------------------- call BeginSub( subname ) if (.not. axis_x_initialized) then call EndSub( subname, 'Call axis_x_init before call %c', & & c1=trim(subname) ) return endif !---------------------------------------------------------------- ! decision が '***@lon' などでない場合は停止して返す。 !---------------------------------------------------------------- if (axis_x_data_from_netcdf == '') then call EndSub( subname, & & 'In axis_x_nml, configurated to not generate from netcdf data') return endif !---------------------------------------------------------------- ! HistoryGet によるデータ入力 !---------------------------------------------------------------- call UrlSplit(axis_x_data_from_netcdf, file=file, var=varname) call HistoryGet(file, varname, x_Dim%a_Dim) x_Dim%stored = .true. ! 入力するデータの units が 'degree' か 'radian' か調べ、 ! 一時、'radian' に変換し、その後、x_Dim%axisinfo%units が ! 'degree' か 'radian' かに応じて 180. / pi を掛けるかどうか ! 判定すべき !---------------------------------------------------------------- ! intent(out) の引数にデータを渡す。 !---------------------------------------------------------------- call axis_type_copy(x_Dim, Dim) call EndSub( subname ) end subroutine
スペクトル法を用いる場合を想定した X 軸のデータを返す。
((< axis_x_init >)) の NAMELIST axis_x_nml の units に (({ ‘radian’ })) または (({ ‘rad.’ })) を与える場合には 単位がラジアンでデータが返される。それ以外では度数でデータが返る。
((< axis_x_init >)) の NAMELIST axis_x_nml の decision 変数で (({ ‘spectral’ })) 以外が与えられた場合は値を代入しないで返す。
subroutine axis_x_spectral (Dim) !==== Dependency use axis_type_mod, only: axis_type_copy use constants_mod, only: constants_init, pi use spml_mod, only: wa_module_x_Lon => x_Lon use dc_string, only: toChar, StrHead, LChar use dc_trace, only: DbgMessage, BeginSub, EndSub !=end implicit none !=begin !==== In/Out ! type(AXISINFO), intent(inout) :: Dim ! 次元情報を包括する変数 !=end character(STRING), parameter:: subname = "axis_x_spectral" continue !---------------------------------------------------------------- ! Check Initialization !---------------------------------------------------------------- call BeginSub( subname ) if (.not. axis_x_initialized) then call EndSub( subname, 'Call axis_x_init before call %c', & & c1=trim(subname) ) return endif !---------------------------------------------------------------- ! decision が 'spectral' でない場合は停止して返す。 !---------------------------------------------------------------- if (.not. axis_x_data_from_spectral) then call EndSub( subname, & & 'In axis_x_nml, configurated to not generate from spmodel') return endif !---------------------------------------------------------------- ! Initialize Dependent module !---------------------------------------------------------------- call constants_init() !---------------------------------------------------------------- ! wa_module からスペクトルデータを受け取る。 !---------------------------------------------------------------- if ( StrHead( 'radians', trim(LChar(x_Dim%axisinfo%units)) ) .or.& & StrHead( 'rad.', trim(LChar(x_Dim%axisinfo%units)) ) ) then ! radian でそのまま代入 x_Dim%a_Dim(:) = wa_module_x_Lon(:) else ! radian を degree に変換して代入 x_Dim%a_Dim(:) = wa_module_x_Lon(:) * 180. / pi end if x_Dim%stored = .true. call DbgMessage('x_Lon=<%c>', c1=trim(toChar(x_Dim%a_Dim)) ) !---------------------------------------------------------------- ! intent(out) の引数にデータを渡す。 !---------------------------------------------------------------- call axis_type_copy(x_Dim, Dim) call EndSub( subname ) end subroutine
重みデータとその付加情報を返す。 ((< axis_x_init >)) の NAMELIST axis_x_nml の decision 変数で (({ ‘spectral’ })) 以外が与えられた場合は値を代入しないで返す。
また、X 軸の次元変数に重みデータに関する付加情報を加える。
subroutine axis_x_weight (Dim_Weight) !==== Dependency use constants_mod, only: constants_init, pi use axis_type_mod, only: axis_type_copy, axis_attrs_copy, axis_attrs_init use spml_mod, only: wa_module_x_Lon_Weight => x_Lon_Weight use grid_3d_mod,only: im use gt4_history,only: GT_HISTORY_ATTR use dc_trace, only: DbgMessage, BeginSub, EndSub !=end implicit none !=begin !==== Output ! type(AXISINFO), intent(out) :: Dim_Weight ! 次元情報を包括する変数 !=end type(GT_HISTORY_ATTR), allocatable :: attrs_tmp(:) ! 次元情報を一時格納 character(STRING), parameter:: subname = "axis_x_weight" continue !---------------------------------------------------------------- ! Check Initialization !---------------------------------------------------------------- call BeginSub( subname ) if (.not. axis_x_initialized) then call EndSub( subname, 'Call axis_x_init before call %c', & & c1=trim(subname) ) return endif !---------------------------------------------------------------- ! decision が 'spectral' でない場合は停止して返す。 !---------------------------------------------------------------- if (.not. axis_x_data_from_spectral) then call EndSub( subname, & & 'In axis_x_nml, configurated to not generate from spmodel') return endif !---------------------------------------------------------------- ! Initialize Dependent module !---------------------------------------------------------------- call constants_init() !---------------------------------------------------------------- ! x_Dim の情報から x_Dim_Weight の情報を生成。 !---------------------------------------------------------------- x_Dim_Weight%axisinfo%name = trim(x_Dim%axisinfo%name) // '_weight' x_Dim_Weight%axisinfo%length = x_Dim%axisinfo%length x_Dim_Weight%axisinfo%longname = trim(x_Dim%axisinfo%longname) // ' weight' x_Dim_Weight%axisinfo%units = trim(x_Dim%axisinfo%units) x_Dim_Weight%axisinfo%xtype = trim(x_Dim%axisinfo%xtype) if ( allocated(x_Dim_Weight%attrs) ) then deallocate( x_Dim_Weight%attrs ) endif ! x_Dim にも情報付加 if ( allocated(x_Dim%attrs) ) then allocate( attrs_tmp(size(x_Dim%attrs)) ) call axis_attrs_copy(from=x_Dim%attrs(:), to=attrs_tmp(:)) deallocate( x_Dim%attrs ) allocate( x_Dim%attrs(size(attrs_tmp)+1) ) call axis_attrs_copy(from=attrs_tmp(:), to=x_Dim%attrs(1:size(attrs_tmp)) ) call axis_attrs_init( x_Dim%attrs(size(attrs_tmp)+1) ) x_Dim%attrs(size(attrs_tmp)+1)%attrname = 'gt_calc_weight' x_Dim%attrs(size(attrs_tmp)+1)%attrtype = 'char' x_Dim%attrs(size(attrs_tmp)+1)%cvalue = x_Dim_Weight%axisinfo%name x_Dim%attrs(size(attrs_tmp)+1)%array = .false. else allocate( x_Dim%attrs(1) ) call axis_attrs_init( x_Dim%attrs(1) ) x_Dim%attrs(1)%attrname = 'gt_calc_weight' x_Dim%attrs(1)%attrtype = 'char' x_Dim%attrs(1)%cvalue = x_Dim_Weight%axisinfo%name x_Dim%attrs(1)%array = .false. endif !---------------------------------------------------------------- ! wa_module からスペクトルデータを受け取る。 !---------------------------------------------------------------- allocate( x_Dim_Weight%a_Dim(im) ) ! ラジアンを度数に変換して代入 x_Dim_Weight%a_Dim(:) = wa_module_x_Lon_Weight(:) * 180. / pi x_Dim_Weight%stored = .true. !---------------------------------------------------------------- ! intent(out) の引数にデータを渡す。 !---------------------------------------------------------------- call axis_type_copy(x_Dim_Weight, Dim_Weight) call EndSub( subname ) end subroutine