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
! 格子点数・最大波数設定
! 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
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
#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, OMPNumThreads, 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
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( r_SSDepth( 0:kslmax ) )
allocate( z_SSDepth( 1:kslmax ) )
! デフォルト値の設定
! Default values settings
!
OMPNumThreads = 1
!$ OMPNumThreads = omp_get_max_threads()
! 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)
!
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)
! 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 ) ! (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
! スペクトルデータの添字番号の設定
! 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