| Class | initial_data |
| In: |
prepare_data/initial_data.f90
|
初期値データ (リスタートデータ) のサンプルを提供します.
現在は以下のデータを提供します.
Prepare sample data of initial data (restart data)
Now, following data are provided.
| InitDataGet : | 初期値データ (リスタートデータ) 取得 |
| ———— : | ———— |
| InitDataGet : | Get initial data (restart data) |
| Subroutine : | |||
| xyz_UB(0:imax-1, 1:jmax, 1:kmax) : | real(DP), intent(out)
| ||
| xyz_VB(0:imax-1, 1:jmax, 1:kmax) : | real(DP), intent(out)
| ||
| xyz_TempB(0:imax-1, 1:jmax, 1:kmax) : | real(DP), intent(out)
| ||
| xyz_QVapB(0:imax-1, 1:jmax, 1:kmax) : | real(DP), intent(out)
| ||
| xy_PsB(0:imax-1, 1:jmax) : | real(DP), intent(out)
| ||
| xyz_UN(0:imax-1, 1:jmax, 1:kmax) : | real(DP), intent(out)
| ||
| xyz_VN(0:imax-1, 1:jmax, 1:kmax) : | real(DP), intent(out)
| ||
| xyz_TempN(0:imax-1, 1:jmax, 1:kmax) : | real(DP), intent(out)
| ||
| xyz_QVapN(0:imax-1, 1:jmax, 1:kmax) : | real(DP), intent(out)
| ||
| xy_PsN(0:imax-1, 1:jmax) : | real(DP), intent(out)
|
初期値データ (リスタートデータ) のサンプルを提供します.
Prepare sample data of initial data (restart data)
subroutine InitDataGet( xyz_UB, xyz_VB, xyz_TempB, xyz_QVapB, xy_PsB, xyz_UN, xyz_VN, xyz_TempN, xyz_QVapN, xy_PsN )
!
! 初期値データ (リスタートデータ) のサンプルを提供します.
!
! Prepare sample data of initial data (restart data)
!
! モジュール引用 ; USE statements
!
! 座標データ設定
! Axes data settings
!
use axesset, only: x_Lon, y_Lat, z_Sigma
! $ \sigma $ レベル (整数).
! Full $ \sigma $ level
! 物理定数設定
! Physical constants settings
!
use constants, only: CpDry, GasRDry ! $ R $ [J kg-1 K-1].
! 乾燥大気の気体定数.
! Gas constant of air
! 飽和比湿計算 (Nakajima et al., 1992)
! Evaluate saturation specific humidity (Nakajima et al., 1992)
!
use saturate_nha1992, only: CalcQVapSat
! 文字列操作
! Character handling
!
use dc_string, only: LChar
! 宣言文 ; Declaration statements
!
implicit none
real(DP), intent(out):: xyz_UB (0:imax-1, 1:jmax, 1:kmax)
! $ u (t-\Delta t) $ . 東西風速. Eastward wind
real(DP), intent(out):: xyz_VB (0:imax-1, 1:jmax, 1:kmax)
! $ v (t-\Delta t) $ . 南北風速. Northward wind
real(DP), intent(out):: xyz_TempB (0:imax-1, 1:jmax, 1:kmax)
! $ T (t-\Delta t) $ . 温度. Temperature
real(DP), intent(out):: xyz_QVapB (0:imax-1, 1:jmax, 1:kmax)
! $ q (t-\Delta t) $ . 比湿. Specific humidity
real(DP), intent(out):: xy_PsB (0:imax-1, 1:jmax)
! $ p_s (t-\Delta t) $ . 地表面気圧. Surface pressure
real(DP), intent(out):: xyz_UN (0:imax-1, 1:jmax, 1:kmax)
! $ u (t) $ . 東西風速. Eastward wind
real(DP), intent(out):: xyz_VN (0:imax-1, 1:jmax, 1:kmax)
! $ v (t) $ . 南北風速. Northward wind
real(DP), intent(out):: xyz_TempN (0:imax-1, 1:jmax, 1:kmax)
! $ T (t) $ . 温度. Temperature
real(DP), intent(out):: xyz_QVapN (0:imax-1, 1:jmax, 1:kmax)
! $ q (t) $ . 比湿. Specific humidity
real(DP), intent(out):: xy_PsN (0:imax-1, 1:jmax)
! $ p_s (t) $ . 地表面気圧. Surface pressure
! Sugiyama et al. (2008) 用作業変数
! Work variables for Sugiyama et al. (2008)
!
real(DP):: xyz_PotTemp (0:imax-1, 1:jmax, 1:kmax)
! 温位. Potential temperature
real(DP):: xyz_Press (0:imax-1, 1:jmax, 1:kmax)
! 気圧. Air pressure
real(DP):: xy_TempMin (0:imax-1, 1:jmax)
! 温度の最小値. Minimum value of temperature
real(DP):: xyz_QVapSat (0:imax-1, 1:jmax, 1:kmax)
! 飽和比湿. Saturation specific humidity
! 作業変数
! Work variables
!
integer:: i ! 経度方向に回る DO ループ用作業変数
! Work variables for DO loop in longitude
integer:: j ! 緯度方向に回る DO ループ用作業変数
! Work variables for DO loop in latitude
integer:: k ! 鉛直方向に回る DO ループ用作業変数
! Work variables for DO loop in vertical direction
! 実行文 ; Executable statement
if ( .not. initial_data_inited ) call InitDataInit
! 微小な温度擾乱のある静止場
! Stationary field with small disturbance of temperature
!
select case ( LChar( trim(Pattern) ) )
case ( 'small disturbance of temperature' )
xyz_UB = 0.0_DP
xyz_VB = 0.0_DP
xyz_TempB = TempAvr
xy_PsB = PsAvr
xyz_QVapB = QVapAvr
! 温度に擾乱を与える
! Add perturbation to temperature
!
do k = 1, kmax
do j = 1, jmax
do i = 0, imax - 1
xyz_TempB(i,j,k) = xyz_TempB(i,j,k) + 0.1_DP * sin ( x_Lon(i) * y_Lat(j) ) - 0.1_DP * ( 1.0_DP - z_Sigma(k) )
end do
end do
end do
! ステップ $ t-\Delta t $ のデータをステップ $ t $ へコピー
! Copy data on step $ t-\Delta t $ to step $ t $
!
xyz_UN = xyz_UB
xyz_VN = xyz_VB
xyz_TempN = xyz_TempB
xy_PsN = xy_PsB
xyz_QVapN = xyz_QVapB
! AGCM5.3 のデフォルト初期値
! AGCM5.3 default initial values
!
case ( 'agcm 5.3 default' )
xyz_UB = 0.0_DP
xyz_VB = 0.0_DP
xyz_TempB = TempAvr
xy_PsB = PsAvr
xyz_QVapB = QVapAvr
! 温度に擾乱を与える
! Add perturbation to temperature
!
do k = 1, kmax
do j = 1, jmax
do i = 0, imax - 1
xyz_TempB(i,j,k) = xyz_TempB(i,j,k) + 0.1_DP * sin ( real( ( i + 1 ) * ( jmax - j + 1 ) * ( kmax - k ), DP ) / real( imax * jmax * kmax, DP ) * 10.0_DP )
end do
end do
end do
! ステップ $ t-\Delta t $ のデータをステップ $ t $ へコピー
! Copy data on step $ t-\Delta t $ to step $ t $
!
xyz_UN = xyz_UB
xyz_VN = xyz_VB
xyz_TempN = xyz_TempB
xy_PsN = xy_PsB
xyz_QVapN = xyz_QVapB
! Sugiyama et al. (2008) の初期値
! Initial values of Sugiyama et al. (2008)
!
case ( 'sugiyama et al. (2008)' )
xyz_UB = 0.0_DP
xyz_VB = 0.0_DP
xyz_TempB = TempAvr
xy_PsB = PsAvr
xyz_QVapB = QVapAvr
! 温度の計算
! Calculate temperature
!
xyz_PotTemp = TempAvr
xy_TempMin = TempAvr
do k = 1, kmax
xyz_TempB(:,:,k) = xyz_PotTemp(:,:,k) * ( z_Sigma(k) )**( GasRDry / CpDry )
if ( PsAvr * z_Sigma(k) < 1.0e+4_DP ) then
xyz_TempB(:,:,k) = xy_TempMin
else
xy_TempMin = xyz_TempB(:,:,k)
end if
end do
! 温度に擾乱を与える
! Add perturbation to temperature
!
do k = 1, kmax
do j = 1, jmax
do i = 0, imax - 1
xyz_TempB(i,j,k) = xyz_TempB(i,j,k) + 0.1_DP * sin ( real( ( i + 1 ) * ( jmax - j + 1 ) * ( kmax - k ), DP ) / real( imax * jmax * kmax, DP ) * 10.0_DP )
end do
end do
end do
! 飽和比湿計算
! Calculate saturation specific humidity
!
do k = 1, kmax
xyz_Press(:,:,k) = xy_PsB * z_Sigma(k)
end do
call CalcQVapSat( xyz_TempB, xyz_Press, xyz_QVapSat ) ! (out)
!!$ call PhySatNhaCreate( &
!!$ & phy_sat_nha = phy_sat, & ! (out)
!!$ & imax = imax, jmax = jmax, kmax = kmax, & ! (in)
!!$ & EpsV = EpsV, & ! (in)
!!$ & err = err ) ! (out)
!!$
!!$ do k = 0, kmax - 1
!!$ xyz_Press(:,:,k) = xy_Ps * z_Sigma(k)
!!$ end do
!!$
!!$ call CalcQVapSat( &
!!$ & phy_sat_nha = phy_sat, & ! (in)
!!$ & xyz_Temp = xyz_Temp, & ! (in)
!!$ & xyz_Press = xyz_Press, & ! (in)
!!$ & xyz_QVapSat = xyz_QVapSat, & ! (out)
!!$ & err = err ) ! (out)
! 比湿の計算
! Calculate specific humidity
!
where ( xyz_QVapB > xyz_QVapSat * 0.75 )
xyz_QVapB = xyz_QVapSat * 0.75
end where
! ステップ $ t-\Delta t $ のデータをステップ $ t $ へコピー
! Copy data on step $ t-\Delta t $ to step $ t $
!
xyz_UN = xyz_UB
xyz_VN = xyz_VB
xyz_TempN = xyz_TempB
xy_PsN = xy_PsB
xyz_QVapN = xyz_QVapB
end select
end subroutine InitDataGet
| Variable : | |||
| initial_data_inited = .false. : | logical, save, public
|
| Subroutine : |
依存モジュールの初期化チェック
Check initialization of dependency modules
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
! 座標データ設定
! Axes data settings
!
use axesset, only: axesset_inited
! 実行文 ; 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.' )
if ( .not. axesset_inited ) call MessageNotify( 'E', module_name, '"axesset" module is not initialized.' )
end subroutine InitCheck
| Subroutine : |
This procedure input/output NAMELIST#initial_data_nml .
subroutine InitDataInit
! モジュール引用 ; USE statements
!
! NAMELIST ファイル入力に関するユーティリティ
! Utilities for NAMELIST file input
!
use namelist_util, only: namelist_filename, NmlutilMsg, NmlutilAryValid
! ファイル入出力補助
! File I/O support
!
use dc_iounit, only: FileOpen
! 文字列操作
! Character handling
!
use dc_string, only: LChar
! 宣言文 ; Declaration statements
!
implicit none
integer:: unit_nml ! NAMELIST ファイルオープン用装置番号.
! Unit number for NAMELIST file open
integer:: iostat_nml ! NAMELIST 読み込み時の IOSTAT.
! IOSTAT of NAMELIST read
! NAMELIST 変数群
! NAMELIST group name
!
namelist /initial_data_nml/ Pattern, TempAvr, PsAvr, QVapAvr
!
! デフォルト値については初期化手続 "initial_data#InitDataInit"
! のソースコードを参照のこと.
!
! Refer to source codes in the initialization procedure
! "initial_data#InitDataInit" for the default values.
!
! 実行文 ; Executable statement
if ( initial_data_inited ) return
call InitCheck
! デフォルト値の設定 (まずは Pattern のみ)
! Default values settings (At first, "Pattern" only)
!
Pattern = 'Small Disturbance of Temperature'
!Pattern = 'AGCM 5.3 Default'
! NAMELIST の読み込み (まずは Pattern のみ)
! NAMELIST is input (At first, "Pattern" only)
!
if ( trim(namelist_filename) /= '' ) then
call FileOpen( unit_nml, namelist_filename, mode = 'r' ) ! (in)
rewind( unit_nml )
read( unit_nml, nml = initial_data_nml, iostat = iostat_nml ) ! (out)
close( unit_nml )
call NmlutilMsg( iostat_nml, module_name ) ! (in)
end if
! デフォルト値の設定
! Default values settings
!
select case ( LChar( trim(Pattern) ) )
case ( 'small disturbance of temperature' )
TempAvr = 250.0_DP
PsAvr = 1.0e+5_DP
QVapAvr = 1.0e-10_DP
case ( 'agcm 5.3 default' )
TempAvr = 250.0_DP
PsAvr = 1.0e+5_DP
QVapAvr = 1.0e-10_DP
case ( 'sugiyama et al. (2008)' )
TempAvr = 490.0_DP
PsAvr = 3.0e+6_DP
QVapAvr = 6.11641e-3_DP
case default
call MessageNotify( 'E', module_name, 'Pattern=<%c> is invalid.', c1 = trim(Pattern) )
end select
! 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 = initial_data_nml, iostat = iostat_nml ) ! (out)
close( unit_nml )
call NmlutilMsg( iostat_nml, module_name ) ! (in)
end if
! 印字 ; Print
!
call MessageNotify( 'M', module_name, '----- Initialization Messages -----' )
call MessageNotify( 'M', module_name, ' Pattern = %c', c1 = trim(Pattern) )
call MessageNotify( 'M', module_name, ' TempAvr = %f', d = (/ TempAvr /) )
call MessageNotify( 'M', module_name, ' PsAvr = %f', d = (/ PsAvr /) )
call MessageNotify( 'M', module_name, ' QVapAvr = %f', d = (/ QVapAvr /) )
call MessageNotify( 'M', module_name, '-- version = %c', c1 = trim(version) )
initial_data_inited = .true.
end subroutine InitDataInit
| Constant : | |||
| module_name = ‘initial_data‘ : | character(*), parameter
|
| Constant : | |||
| version = ’$Name: dcpam5-20081109-1 $’ // ’$Id: initial_data.f90,v 1.4 2008-10-21 08:55:58 morikawa Exp $’ : | character(*), parameter
|