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 :経度座標重み
DeltaLon :経度格子点間隔
InvDeltaLon :経度格子点間隔の逆数
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
DeltaLon :Interval of longitude grids
InvDeltaLon :Inverse of DeltaLon
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

Procedures List

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

Methods

Included Modules

gridset dc_types namelist_util constants0 constants wa_mpi_module_sjpack wa_mpi_module wa_zonal_module wa_module_sjpack wa_zonal_module_sjpack 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. 
    !

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

    ! 宣言文 ; 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( r_SSDepth ) ) deallocate( r_SSDepth )
    if ( allocated( z_SSDepth ) ) deallocate( z_SSDepth )

    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 and mathematical constants settings
    !
    use constants0, only: PI      ! $ \pi $ .
                                  ! 円周率.  Circular constant

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

    ! 格子点数・最大波数設定
    ! Number of grid points and maximum truncated wavenumber settings
    !
    use gridset, only: GridsetCheckNumberOfLatGrid

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

    ! OpenMP
    !
    !$ use omp_lib


    ! 宣言文 ; Declaration statements
    !
    implicit none

    ! 作業変数
    ! Work variables
    !
    logical:: flag_generate_sigma
                              ! 鉛直層数の内部生成のためのフラグ
                              ! Flag for generation of sigma levels internally
    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

    real(DP):: LonInDeg       ! Longitude in unit of degree of a grid point in 1D mode
    real(DP):: LatInDeg       ! Latitude in unit of degree of a grid point in 1D mode

    integer:: unit_nml        ! NAMELIST ファイルオープン用装置番号. 
                              ! Unit number for NAMELIST file open
    integer:: iostat_nml      ! NAMELIST 読み込み時の IOSTAT. 
                              ! IOSTAT of NAMELIST read

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

    character(TOKEN):: rank_str
                              ! ランク数. Rank number
    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, Depth, flag_generate_sigma, LonInDeg, LatInDeg
          !
          ! デフォルト値については初期化手続 "axesset#AxessetInit" 
          ! のソースコードを参照のこと. 
          !
          ! Refer to source codes in the initialization procedure
          ! "axesset#AxessetInit" for the default values. 
          !

    ! 実行文 ; Executable statement
    !

    if ( axesset_inited ) return


    !
    ! Set number of OpenMP threads
    !
    OMPNumThreads = 1
    !$ OMPNumThreads = omp_get_max_threads()


    ! 割り付け
    ! 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( r_SSDepth( 0:kslmax ) )
    allocate( z_SSDepth( 1:kslmax ) )

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

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

    ! 地下の層の境界面深さの初期値 (無効な値) の設定
    ! Setting of initial value (invalid value) of depth of subsurface layer interface
    !
    Depth = -999.0d0

    ! 鉛直層数の内部生成のためのフラグの設定
    ! Setting of flag for generation of sigma levels internally
    !
    flag_generate_sigma = .false.

    ! Longitude in unit of degree of a grid point in 1D mode
    !
    LonInDeg = 0.0_DP

    ! Latitude in unit of degree of a grid point in 1D mode
    !
    LatInDeg  = 0.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 == -999.0d0 ) ) then
      if ( flag_generate_sigma ) then
        call SigmaDataGetHalf( Sigma(1:kmax+1) ) ! (out)
      else
        call MessageNotify( 'E', module_name, '  Sigma levels have to be specified if flag_generate_sigma is not true.' )
      end if
    end if

    ! Sigma (半整数レベルσ) チェック
    ! Check "Sigma" (half level sigma)
    !
    if ( Sigma(1) /= 1.0_DP ) then
      call MessageNotify( 'E', module_name, '  Sigma(1) is not 1, but is %f', d = (/ Sigma(1) /) )
    end if
    if ( Sigma(kmax+1) /= 0.0_DP ) then
      call MessageNotify( 'E', module_name, '  Sigma(kmax+1) is not 0, but is %f', d = (/ Sigma(kmax+1) /) )
    end if
    do k = 1, kmax
      if ( Sigma(k+1) > Sigma(k) ) then
        call MessageNotify( 'M', module_name, '  Sigma(%d) = %f > Sigma(%d) = %f', i = (/ k+1, k /), d = (/ Sigma(k+1), Sigma(k) /) )
        call MessageNotify( 'E', module_name, '  Value of Sigma has to decrease with index.' )
      end if
    end do

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

    ! 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

    ! 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


    ! 地下の層の境界面深さのチェック
    ! Check depth of subsurface layer interface
    !
    if ( all( Depth == -999.0d0 ) ) then
      Depth(0+1:kslmax+1 ) = 0.0d0
      call MessageNotify( 'W', module_name, 'Depth is not found in namelist file.' )
    end if
    !
    if ( Depth(0+1) /= 0.0d0 ) then
      call MessageNotify( 'E', module_name, '  Depth(0) is not zero, but is %f', d = (/ Depth(0+1) /) )
    end if
    if ( kslmax >= 1 ) then
      if ( all( Depth(1+1:kslmax+1) >= 0.0d0 ) ) then
        do k = 0, kslmax
          call MessageNotify( 'M', module_name, '  Depth(%d) = %f', i = (/ k /), d = (/ Depth(k+1) /) )
        end do
        call MessageNotify( 'E', module_name, '  Depth has to be zero or negative.' )
      end if
    end if
    do k = 0, kslmax-1
      if ( Depth(k+1+1) > Depth(k+1) ) then
        call MessageNotify( 'M', module_name, '  Depth(%d) = %f > Depth(%d) = %f', i = (/ k+1, k /), d = (/ Depth(k+1+1), Depth(k+1) /) )
        call MessageNotify( 'E', module_name, '  Value of Depth has to decrease with index.' )
      end if
    end do

    r_SSDepth(0:kslmax) = Depth(1:kslmax+1)

    do k = 0, kslmax
      call MessageNotify( 'M', module_name, '  r_SSDepth(%d) = %f', i = (/ k /), d = (/ r_SSDepth(k) /) )
    end do

    ! 地下の鉛直層の中心点の設定
    ! Set midpoint of subsurface grid
    !
    do k = 1, kslmax
      z_SSDepth( k ) = ( r_SSDepth( k-1 ) + r_SSDepth( k ) ) / 2.0d0
    end do

    ! 緯度経度の設定
    ! Settings of longitude and latitude
    !
    allocate( x_Lon        (0:imax-1) )
    allocate( x_Lon_Weight (0:imax-1) )
    allocate( y_Lat        (1:jmax)   )
    allocate( y_Lat_Weight (1:jmax)   )

    if ( ( imax == 1 ) .and. ( jmax == 1 ) ) then

      x_Lon          = LonInDeg * PI / 180.0_DP
      x_Lon_Weight   = 1.0_DP
      y_Lat          = LatInDeg * PI / 180.0_DP
      y_Lat_Weight   = 1.0_DP

    else
      if ( .not. spml_inited ) then
#ifdef LIB_MPI
        call wa_mpi_Initial( nmax, imax, jmax_global, kmax, OMPNumThreads ) ! (in)
        ! Check number of latitudinal grid on each process
        call GridsetCheckNumberOfLatGrid( spml_jc )
#else
        call wa_Initial( nmax, imax, jmax_global, kmax, OMPNumThreads ) ! (in)
#endif
        spml_inited = .true.
      end if

      x_Lon          = spml_x_Lon
      x_Lon_Weight   = spml_x_Lon_Weight
      y_Lat          = spml_y_Lat
      y_Lat_Weight   = spml_y_Lat_Weight
    end if

    if ( imax /= 1  ) then
      ! DeltaLon, InvDeltaLonの計算
      ! Caluculate DeltaLon and InvDeltaLon
      DeltaLon = x_Lon(1) - x_Lon(0)
      InvDeltaLon = 1.0_DP/DeltaLon
    else
      DeltaLon = 0.0_DP
      InvDeltaLon = 0.0_DP ! not used
    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
    rank_str = ''


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

#ifdef SJPACK
    call MessageNotify( 'M', module_name, 'SJPACK : %c', c1 = 'used.' )
#else
    call MessageNotify( 'M', module_name, 'SJPACK : %c', c1 = 'not used.' )
#endif

    call MessageNotify( 'M', module_name, 'OMPNumThreads = %d', i = (/OMPNumThreads/) )

    do ra = 0, nprocs_mpi - 1

      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
AxnameR
Constant :
AxnameR = ‘sigm‘ :character(STRING), parameter, public
: axis name for r
AxnameSSR
Constant :
AxnameSSR = ‘sszi‘ :character(STRING), parameter, public
: axis name for r (subsurface)
AxnameSSZ
Constant :
AxnameSSZ = ‘ssz‘ :character(STRING), parameter, public
: axis name for z (subsurface)
AxnameT
Constant :
AxnameT = ‘time‘ :character(STRING), parameter, public
: axis name for t
AxnameWN
Constant :
AxnameWN = ‘wn‘ :character(STRING), parameter, public
: axis name for wavenumber
AxnameX
Constant :
AxnameX = ‘lon‘ :character(STRING), parameter, public
: axis name for x
AxnameY
Constant :
AxnameY = ‘lat‘ :character(STRING), parameter, public
: axis name for y
AxnameZ
Constant :
AxnameZ = ‘sig‘ :character(STRING), parameter, public
: axis name for z
DeltaLon
Variable :
DeltaLon :real(DP), save, public
: 経度格子点間隔. grid interval in longitude
InvDeltaLon
Variable :
InvDeltaLon :real(DP), save, public
: 経度格子点間隔の逆数. Inverse of the grid interval in longitude
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_SSDepth
Variable :
r_SSDepth(:) :real(DP), allocatable, save, public
: 地下の層の境界面の深さ subsurface grid on interface of layer
r_Sigma
Variable :
r_Sigma(:) :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
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
z_DelSigma
Variable :
z_DelSigma(:) :real(DP), allocatable, save, public
: $ Delta sigma $ (整数). $ Delta sigma $ (Full)
z_SSDepth
Variable :
z_SSDepth(:) :real(DP), allocatable, save, public
: 地下の格子点の深さ subsurface grid at midpoint of layer
z_Sigma
Variable :
z_Sigma(:) :real(DP), allocatable, save, public
: $ sigma $ レベル (整数). Full $ sigma $ level

Private Instance methods

Depth
Variable :
Depth(1:MaxNmlArySize) :real(DP), save
: 入力用配列, 地下の鉛直層境界 array for input, layer interface of subsurface vertical layer
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-20140314 $’ // ’$Id: axesset.F90,v 1.19 2013-09-16 12:09:34 yot Exp $’ :character(*), parameter
: モジュールのバージョン Module version