Class axesset
In: setup/axesset.F90

座標データ設定

Axes data settings

Note that Japanese and English are described in parallel.

座標データの設定および保管を行います.

緯度 $ varphi $ および経度 $ lambda $ の格子点は, gridset モジュールで設定される格子点数から, SPMODEL ライブラリ を用いて決定されます. 緯度の格子点の位置はガウス緯度, 経度の格子点の位置は等間隔にとることになります.

鉛直σの格子点は, 半整数レベル (各層の端) について NAMELIST#Sigma に指定します. 整数レベル (各層の中心点) は, 半整数レベルと constants モジュールで設定される乾燥大気の定圧比熱 $ C_p $ および 乾燥大気の気体定数 $ R $ を用いて決定します. 鉛直σの格子点については, sigma_data モジュールで用意されている ものを使用することも可能です.

Axes data is set and stored.

Grid points of latitude $ varphi $ and longitude $ lambda $ are determined with number of grid points configured in "gridset" module by SPMODEL Library Grid points of latitude becomes Gaussian latitude, and grid points of latitude becomes equally spaced.

Variables List

x_Lon :経度座標
x_Lon_Weight :経度座標重み
y_Lat :緯度座標
y_Lat_Weight :緯度座標重み
z_Sigma :$ sigma $ レベル (整数)
r_Sigma :$ sigma $ レベル (半整数)
z_DelSigma :$ Delta sigma $ (整数)
r_DelSigma :$ Delta sigma $ (半整数)
w_Number :スペクトルデータの添字番号
spml_inited :SPML ライブラリの初期設定フラグ
———— :————
x_Lon :Longitude
x_Lon_Weight :Weight of longitude
y_Lat :Latitude
y_Lat_Weight :Weight of latitude
z_Sigma :Full $ sigma $ level
r_Sigma :Half $ sigma $ level
z_DelSigma :$ Delta sigma $ (Full)
r_DelSigma :$ Delta sigma $ (Half)
w_Number :Subscript of spectral data
spml_inited :Initialization flag of SPML library

Variables List (MPI)

以下は MPI 使用時に, 領域全体の座標データを参照する際に使用してください.

Use following variables in order to refer data of axes in whole area.

x_Lon_wholeMPI :経度座標
y_Lat_wholeMPI :緯度座標
z_Sigma_wholeMPI :$ sigma $ レベル (整数)
r_Sigma_wholeMPI :$ sigma $ レベル (半整数)
———— :————
x_Lon_wholeMPI :Longitude
y_Lat_wholeMPI :Latitude
z_Sigma_wholeMPI :Full $ sigma $ level
r_Sigma_wholeMPI :Half $ sigma $ level

Procedures List

AxessetInit :座標データの設定
AxessetFinalize :終了処理 (モジュール内部の変数の割り付け解除)
———— :————
AxessetInit :Settings of axes data
AxessetFinalize :Termination (deallocate variables in this module)

Methods

Included Modules

gridset dc_types namelist_util constants wa_mpi_module wa_module sigma_data dc_iounit dc_message dc_string

Public Instance methods

Subroutine :

モジュール内部の変数の割り付け解除を行います.

Deallocate variables in this module.

[Source]

  subroutine AxessetFinalize
    !
    ! モジュール内部の変数の割り付け解除を行います. 
    !
    ! Deallocate variables in this module. 
    !

    ! 宣言文 ; Declaration statements
    !
    implicit none

    ! 実行文 ; Executable statement
    !

    if ( .not. axesset_inited ) return

    ! 割り付け解除
    ! Deallocation
    !
    if ( allocated( z_Sigma    ) ) deallocate( z_Sigma    )
    if ( allocated( r_Sigma    ) ) deallocate( r_Sigma    )
    if ( allocated( z_DelSigma ) ) deallocate( z_DelSigma )
    if ( allocated( r_DelSigma ) ) deallocate( r_DelSigma )
    if ( allocated( w_Number   ) ) deallocate( w_Number   )

    if ( allocated( z_Sigma_wholeMPI ) ) deallocate( z_Sigma_wholeMPI )
    if ( allocated( r_Sigma_wholeMPI ) ) deallocate( r_Sigma_wholeMPI )
    if ( allocated( x_Lon_wholeMPI   ) ) deallocate( x_Lon_wholeMPI   )
    if ( allocated( y_Lat_wholeMPI   ) ) deallocate( y_Lat_wholeMPI   )

    if ( allocated( x_Lon        ) ) deallocate( x_Lon        )
    if ( allocated( x_Lon_Weight ) ) deallocate( x_Lon_Weight )
    if ( allocated( y_Lat        ) ) deallocate( y_Lat        )
    if ( allocated( y_Lat_Weight ) ) deallocate( y_Lat_Weight )

    axesset_inited = .false.
  end subroutine AxessetFinalize
Subroutine :

axesset モジュールの初期化を行います. NAMELIST#axesset_nml の読み込みはこの手続きで行われます.

"axesset" module is initialized. NAMELIST#axesset_nml is loaded in this procedure.

This procedure input/output NAMELIST#axesset_nml .

[Source]

  subroutine AxessetInit
    !
    ! axesset モジュールの初期化を行います. 
    ! NAMELIST#axesset_nml の読み込みはこの手続きで行われます. 
    !
    ! "axesset" module is initialized. 
    ! NAMELIST#axesset_nml is loaded in this procedure. 
    !

    ! モジュール引用 ; USE statements
    !

    ! 物理定数設定
    ! Physical constants settings
    !
    use constants, only: PI, CpDry, GasRDry  ! $ R $ [J kg-1 K-1]. 
                                  ! 乾燥大気の気体定数. 
                                  ! Gas constant of air

#ifdef LIB_MPI
    ! MPI 版 SPMODEL ライブラリ, 球面上の問題を球面調和函数変換により解く(多層対応) 
    ! SPMODEL library MPI version, problems on sphere are solved with spherical harmonics (multi layer is supported)
    !
    use wa_mpi_module, only: wa_mpi_Initial, spml_jc           => jc, spml_x_Lon        => x_Lon, spml_y_Lat        => v_Lat, spml_y_Lat_wholeMPI => y_Lat, spml_x_Lon_Weight => x_Lon_Weight, spml_y_Lat_Weight => v_Lat_Weight
                              ! $ \Delta \varphi $ [rad.] . 
                              ! 緯度座標重み. 
                              ! Weight of latitude
#else
    ! SPMODEL ライブラリ, 球面上の問題を球面調和函数変換により解く(多層対応) 
    ! SPMODEL library, problems on sphere are solved with spherical harmonics (multi layer is supported)
    !
    use wa_module, only: wa_Initial, spml_x_Lon        => x_Lon, spml_y_Lat        => y_Lat, spml_x_Lon_Weight => x_Lon_Weight, spml_y_Lat_Weight => y_Lat_Weight
                              ! $ \Delta \varphi $ [rad.] . 
                              ! 緯度座標重み. 
                              ! Weight of latitude
#endif

    ! 鉛直σレベルデータ準備
    ! Prepare vertical sigma level data
    !
    use sigma_data, only: SigmaDataGetHalf

    ! NAMELIST ファイル入力に関するユーティリティ
    ! Utilities for NAMELIST file input
    !
    use namelist_util, only: namelist_filename, NmlutilMsg, NmlutilAryValid

    ! ファイル入出力補助
    ! File I/O support
    !
    use dc_iounit, only: FileOpen

    ! 種別型パラメタ
    ! Kind type parameter
    !
    use dc_types, only: STDOUT, TOKEN      ! キーワード.   Keywords. 

    ! メッセージ出力
    ! Message output
    !
    use dc_message, only: MessageNotify

    ! 文字列操作
    ! Character handling
    !
    use dc_string, only: CPrintf

    ! 宣言文 ; Declaration statements
    !
    implicit none

    ! 作業変数
    ! Work variables
    !
    integer:: i               ! スペクトルの添字番号で回る DO ループ用作業変数
                              ! Work variables for DO loop in subscript of spectral data
    integer:: k               ! 鉛直方向に回る DO ループ用作業変数
                              ! Work variables for DO loop in vertical direction
    real(DP):: Kappa          ! $ \kappa = R / C_p $ .
                              ! 乾燥大気における, 気体定数の定圧比熱に対する比.
                              ! Ratio of gas constant to specific heat in dry air
    integer:: unit_nml        ! NAMELIST ファイルオープン用装置番号. 
                              ! Unit number for NAMELIST file open
    integer:: iostat_nml      ! NAMELIST 読み込み時の IOSTAT. 
                              ! IOSTAT of NAMELIST read

    integer:: IspackOpenmpThreads
                              ! OPENMP での最大スレッド数.
                              ! openmp_threads に 1 より大きな値を指定すれば 
                              ! ISPACK[http://www.gfd-dennou.org/library/ispack/] 
                              ! の球面調和函数変換 OPENMP 並列計算
                              ! ルーチンが用いられる. 並列計算を実行するには,
                              ! 実行時に環境変数 OMP_NUM_THREADS 
                              ! を IspackOpenmpThreads 以下の数字に設定する
                              ! 等のシステムに応じた準備が必要となる.
                              !
                              ! IspackOpenmpThreads に 1 より大きな値を
                              ! 指定しなければ並列計算ルーチンは呼ばれない.

    character(TOKEN):: rank_str
                              ! ランク数. Rank number
#ifdef LIB_MPI
    logical:: initflag_mpi
    integer:: err_mpi
#endif
    integer:: myrank_mpi      ! 総プロセス数. Number of total processes
    integer:: nprocs_mpi      ! 自身のプロセス. Number of my process
    integer:: ra              ! MPI のランク数方向に回る DO ループ用作業変数
                              ! Work variables for DO loop in rank number of MPI direction

    ! NAMELIST 変数群
    ! NAMELIST group name
    !
    namelist /axesset_nml/ Sigma, IspackOpenmpThreads
          !
          ! デフォルト値については初期化手続 "axesset#AxessetInit" 
          ! のソースコードを参照のこと. 
          !
          ! Refer to source codes in the initialization procedure
          ! "axesset#AxessetInit" for the default values. 
          !

    ! 実行文 ; Executable statement
    !

    if ( axesset_inited ) return
    call InitCheck

    ! 割り付け
    ! Allocation
    !
    allocate( z_Sigma      (1:kmax)   )
    allocate( r_Sigma      (0:kmax)   )
    allocate( z_DelSigma   (1:kmax)   )
    allocate( r_DelSigma   (0:kmax)   )
    allocate( w_Number     (1:(nmax+1)**2) )

    allocate( z_Sigma_wholeMPI (1:kmax)   )
    allocate( r_Sigma_wholeMPI (0:kmax)   )

    ! デフォルト値の設定
    ! Default values settings
    !
    IspackOpenmpThreads = 0

    ! Sigma (半整数レベルσ) の初期値 (無効な値) の設定
    ! Setting of initial value (invalid value) of "Sigma" (half level sigma)
    !
    Sigma = -999.0_DP

    ! NAMELIST の読み込み
    ! NAMELIST is input
    !
    if ( trim(namelist_filename) /= '' ) then
      call FileOpen( unit_nml, namelist_filename, mode = 'r' ) ! (in)

      rewind( unit_nml )
      read( unit_nml, nml = axesset_nml, iostat = iostat_nml ) ! (out)
      close( unit_nml )

      call NmlutilMsg( iostat_nml, module_name ) ! (in)
    end if

    ! Sigma (半整数レベルσ) の自動設定
    ! Automation setting of "Sigma" (half level sigma)
    !
    if ( all( Sigma < 0.0_DP ) ) then
      call SigmaDataGetHalf( Sigma(1:kmax+1) ) ! (out)
    end if

    ! Sigma (半整数レベルσ) チェック
    ! Check "Sigma" (half level sigma)
    !
    call NmlutilAryValid( module_name, Sigma,  'Sigma', kmax+1, 'kmax+1' )               ! (in)

    ! r_Sigma (半整数レベルσ) 設定
    ! Setting of "r_Sigma" (half level sigma)
    !
    r_Sigma(0:kmax) = Sigma(1:kmax+1)
    r_Sigma_wholeMPI(0:kmax) = r_Sigma(0:kmax)

    ! z_DelSigma (整数レベル $ \Delta \sigma $ ) 設定
    ! Setting of "z_DelSigma" (full level $ \Delta \sigma $ )
    !
    do k = 1, kmax
      z_DelSigma(k) = r_Sigma(k-1) - r_Sigma(k)
    enddo

    ! z_Sigma (整数レベルσ) 設定
    ! Setting of "z_Sigma" (full level sigma)
    !
    Kappa = GasRDry / CpDry
    do k = 1, kmax
      z_Sigma(k) = ( (   r_Sigma(k-1) ** ( 1.0_DP + Kappa ) - r_Sigma(k)   ** ( 1.0_DP + Kappa ) ) / ( z_DelSigma(k) * ( 1.0_DP + Kappa ) ) ) ** ( 1.0_DP / Kappa )
    enddo
    z_Sigma_wholeMPI(1:kmax) = z_Sigma(1:kmax)

    ! r_DelSigma (半整数レベル $ \Delta \sigma $ ) 設定
    ! Setting of "r_DelSigma" (half level $ \Delta \sigma $ )
    !
    r_DelSigma(0)    = r_Sigma(0)    - z_Sigma(1)
    r_DelSigma(kmax) = z_Sigma(kmax) - r_Sigma(kmax)
    do k = 1, kmax - 1
      r_DelSigma(k) = z_Sigma(k) - z_Sigma(k+1)
    end do

    ! 緯度経度の設定
    ! Settings of longitude and latitude
    !
    allocate( x_Lon_wholeMPI (0:imax-1) )
    allocate( y_Lat_wholeMPI (1:jmax)   )

    if ( .not. spml_inited ) then
#ifdef LIB_MPI
      call MPI_Initialized(initflag_mpi, err_mpi)
      if ( .not. initflag_mpi ) then
        call MessageNotify( 'E', module_name, 'Initialize MPI by "MPI_Init" before AxessetInit' )
      end if

      call wa_mpi_Initial( nmax, imax, jmax, kmax ) ! (in)
      call GridsetChange( jm = spml_jc ) ! (in)
#else
      if ( IspackOpenmpThreads < 2 ) then
        call wa_Initial( nmax, imax, jmax, kmax ) ! (in)
      else
        call wa_Initial( nmax, imax, jmax, kmax, IspackOpenmpThreads ) ! (in)
      end if
#endif
      spml_inited = .true.
    end if

    allocate( x_Lon        (0:imax-1) )
    allocate( x_Lon_Weight (0:imax-1) )
    allocate( y_Lat        (1:jmax)   )
    allocate( y_Lat_Weight (1:jmax)   )
    x_Lon          = spml_x_Lon
    x_Lon_wholeMPI = spml_x_Lon
    x_Lon_Weight   = spml_x_Lon_Weight
    y_Lat          = spml_y_Lat
    y_Lat_Weight   = spml_y_Lat_Weight

#ifdef LIB_MPI
    y_Lat_wholeMPI = spml_y_Lat_wholeMPI
#else
    y_Lat_wholeMPI = spml_y_Lat
#endif

    ! スペクトルデータの添字番号の設定
    ! Settings of subscript of spectral data
    !
    do i = 1, size(w_Number)
      w_Number(i) = i
    end do

    ! ランクに関する情報の取得
    ! Get information about rank
    !
    myrank_mpi = -1
    nprocs_mpi = 1
#ifdef LIB_MPI
    call MPI_Initialized(initflag_mpi, err_mpi)
    if ( initflag_mpi ) then
      call MPI_Comm_Rank(MPI_COMM_WORLD, myrank_mpi, err_mpi)
      rank_str = CPrintf( ' [rank=%06d]', i = (/ myrank_mpi /) )
      call MPI_Comm_Size(MPI_COMM_WORLD, nprocs_mpi, err_mpi)
    else
      rank_str = ''
    end if
#else
    rank_str = ''
#endif


    ! 印字 ; Print
    !
    call MessageNotify( 'M', module_name, '----- Initialization Messages -----' )

    do ra = 0, nprocs_mpi - 1

#ifdef LIB_MPI
      if ( initflag_mpi ) call MPI_Barrier(MPI_COMM_WORLD, err_mpi)
      if ( myrank_mpi > -1 .and. ra /= myrank_mpi ) cycle
#endif

      call MessageNotify( 'M', module_name, 'Axes:%c', c1 = trim(rank_str), rank_mpi = -1 )
      call MessageNotify( 'M', module_name, '  x_Lon(%d:%d) [deg.] = %*f', i = (/ 0, imax - 1/), d = format_print(x_Lon / PI * 180.0_DP, imax), n =(/ imax /), rank_mpi = -1 )
      call MessageNotify( 'M', module_name, '  y_Lat(%d:%d) [deg.] = %*f', i = (/ 1, jmax/), d = format_print(y_Lat / PI * 180.0_DP, jmax), n =(/ jmax /), rank_mpi = -1 )
      call MessageNotify( 'M', module_name, '  z_Sigma(%d:%d) = %*f', i = (/ 1, kmax /), d = format_print(z_Sigma, kmax), n =(/ kmax /), rank_mpi = -1 )
      call MessageNotify( 'M', module_name, '  r_Sigma(%d:%d) = %*f', i = (/ 0, kmax /), d = format_print(r_Sigma, kmax+1), n =(/ kmax+1 /), rank_mpi = -1 )
      call MessageNotify( 'M', module_name, '  w_Number(%d:%d) = %d .. %d', i = (/ 1, size(w_Number), 1, size(w_Number) /), rank_mpi = -1 )
  !
      call MessageNotify( 'M', module_name, 'Weight:' )
      call MessageNotify( 'M', module_name, '  x_Lon_Weight(%d:%d) = %*f', i = (/ 0, imax - 1/), d = format_print(x_Lon_Weight, imax), n =(/ imax /), rank_mpi = -1 )
      call MessageNotify( 'M', module_name, '  y_Lat_Weight(%d:%d) = %*f', i = (/ 1, jmax/), d = format_print(y_Lat_Weight, jmax), n =(/ jmax /), rank_mpi = -1 )
      call MessageNotify( 'M', module_name, '  z_DelSigma(%d:%d) = %*f', i = (/ 1, kmax /), d = format_print(z_DelSigma, kmax), n =(/ kmax /), rank_mpi = -1 )
      call MessageNotify( 'M', module_name, '' )
    end do

    call MessageNotify( 'M', module_name, '-- version = %c', c1 = trim(version), rank_mpi = -1 )

    axesset_inited = .true.
  end subroutine AxessetInit
axesset_inited
Variable :
axesset_inited = .false. :logical, save, public
: 初期設定フラグ. Initialization flag
r_DelSigma
Variable :
r_DelSigma(:) :real(DP), allocatable, save, public
: $ Delta sigma $ (半整数). $ Delta sigma $ (half)
r_Sigma
Variable :
r_Sigma(:) :real(DP), allocatable, save, public
: $ sigma $ レベル (半整数). Half $ sigma $ level
r_Sigma_wholeMPI
Variable :
r_Sigma_wholeMPI(:) :real(DP), allocatable, save, public
: $ sigma $ レベル (半整数). Half $ sigma $ level
spml_inited
Variable :
spml_inited = .false. :logical, save, public
: SPML ライブラリの初期設定フラグ. Initialization flag of SPML library
w_Number
Variable :
w_Number(:) :integer, allocatable, save, public
: スペクトルデータの添字番号. Subscript of spectral data
x_Lon
Variable :
x_Lon(:) :real(DP), allocatable, save, public
: $ lambda $ [rad.] . 経度. Longitude
x_Lon_Weight
Variable :
x_Lon_Weight(:) :real(DP), allocatable, save, public
: $ Delta lambda $ [rad.] . 経度座標重み. Weight of longitude
x_Lon_wholeMPI
Variable :
x_Lon_wholeMPI(:) :real(DP), allocatable, save, public
: $ lambda $ [rad.] . 経度. Longitude
y_Lat
Variable :
y_Lat(:) :real(DP), allocatable, save, public
: $ varphi $ [rad.] . 緯度. Latitude
y_Lat_Weight
Variable :
y_Lat_Weight(:) :real(DP), allocatable, save, public
: $ Delta varphi $ [rad.] . 緯度座標重み. Weight of latitude
y_Lat_wholeMPI
Variable :
y_Lat_wholeMPI(:) :real(DP), allocatable, save, public
: $ varphi $ [rad.] . 緯度. Latitude
z_DelSigma
Variable :
z_DelSigma(:) :real(DP), allocatable, save, public
: $ Delta sigma $ (整数). $ Delta sigma $ (Full)
z_Sigma
Variable :
z_Sigma(:) :real(DP), allocatable, save, public
: $ sigma $ レベル (整数). Full $ sigma $ level
z_Sigma_wholeMPI
Variable :
z_Sigma_wholeMPI(:) :real(DP), allocatable, save, public
: $ sigma $ レベル (整数). Full $ sigma $ level

Private Instance methods

Subroutine :

依存モジュールの初期化チェック

Check initialization of dependency modules

[Source]

  subroutine InitCheck
    !
    ! 依存モジュールの初期化チェック
    !
    ! Check initialization of dependency modules

    ! モジュール引用 ; USE statements
    !

    ! NAMELIST ファイル入力に関するユーティリティ
    ! Utilities for NAMELIST file input
    !
    use namelist_util, only: namelist_util_inited

    ! 格子点設定
    ! Grid points settings
    !
    use gridset, only: gridset_inited

    ! 物理定数設定
    ! Physical constants settings
    !
    use constants, only: constants_inited

    ! メッセージ出力
    ! Message output
    !
    use dc_message, only: MessageNotify

    ! 実行文 ; Executable statement
    !

    if ( .not. namelist_util_inited ) call MessageNotify( 'E', module_name, '"namelist_util" module is not initialized.' )

    if ( .not. gridset_inited ) call MessageNotify( 'E', module_name, '"gridset" module is not initialized.' )

    if ( .not. constants_inited ) call MessageNotify( 'E', module_name, '"constants" module is not initialized.' )

  end subroutine InitCheck
Sigma
Variable :
Sigma(1:MaxNmlArySize) :real(DP), save
: $ sigma $ レベル (半整数). Half $ sigma $ level
Function :
result(size) :real(DP)
ary(size) :real(DP), intent(in)
size :integer, intent(in)

標準出力して見やすいように適当に有効数字の桁を落として返す.

Effective digit is reduced for output to standard output and returned.

[Source]

  function format_print( ary, size ) result(result)
    !
    ! 標準出力して見やすいように適当に有効数字の桁を落として返す. 
    !
    ! Effective digit is reduced for output to standard output and returned. 
    !

    ! 宣言文 ; Declaration statements
    !
    implicit none
    integer, intent(in):: size
    real(DP), intent(in):: ary(size)
    real(DP):: result(size)
  
    ! 実行文 ; Executable statement
    !
    result(1:size) = int( ary(1:size) * 1000.0_DP ) / 1000.0_DP
    
  end function format_print
module_name
Constant :
module_name = ‘axesset :character(*), parameter
: モジュールの名称. Module name
version
Constant :
version = ’$Name: dcpam5-20090317 $’ // ’$Id: axesset.F90,v 1.4 2009-03-16 05:28:16 morikawa Exp $’ :character(*), parameter
: モジュールのバージョン Module version

[Validate]