Class initial_data
In: prepare_data/initial_data.f90

初期値データ提供

Prepare initial data

初期値データのサンプルを提供します.

現在は以下のデータを提供します.

  • Small Disturbance of Temperature
    • 風速: 0 [m/s], 地表面気圧: 1.0e+5 [Pa], 比湿: 1.0e-10 [kg kg-1]
    • 温度: 250 [K] に微小擾乱を加えたもの
  • [AGCM 5.3]{www.gfd-dennou.org/library/agcm5} Default
    • 風速: 0 [m/s], 地表面気圧: 1.0e+5 [Pa], 比湿: 1.0e-10 [kg kg-1]
    • 温度: 250 [K] に微小擾乱を加えたもの
  • Sugiyama et al. (2008)
    • 木星大気を想定した, Sugiyama et al. (2008) で用いられていた 初期値を模倣する.
    • 風速: 0 [m/s], 地表面気圧: 3.0e+6 [Pa]
    • 温度: $ sigma = 1 $ で 490 [K] とし, 気圧 1.0e+4 となる高度まで, 温位一定で (乾燥断熱線に沿って) 高度に伴い温度を減少させる. 気圧 1.0e+4 以下の高度では温度一定とする.
    • 比湿: $ sigma = 1 $ で 6.11641e-3 [kg kg-1] とし, 比湿が飽和比湿の 75 % となる高度までは一定とする. その高度以上では, 比湿を飽和比湿の 75 % とする.

Prepare sample data of initial data (restart data)

Now, following data are provided.

  • Small Disturbance of Temperature
    • Velocity: 0 [m/s], Surface pressure: 1.0e+5 [Pa], Specific humidity: 1.0e-10 [kg kg-1]
    • Temperature: 250 [K] and perturbation
  • [AGCM 5.3]{www.gfd-dennou.org/library/agcm5} Default
    • Velocity: 0 [m/s], Surface pressure: 1.0e+5 [Pa], Specific humidity: 1.0e-10 [kg kg-1]
    • Temperature: 250 [K] and perturbation
  • Sugiyama et al. (2008)
    • Initial data like a Jovian atmosphere that is used in Sugiyama et al. (2008) is imitated.
    • Velocity: 0 [m/s], Surface pressure: 3.0e+6 [Pa]
    • Temperature: 490 [K] at $ sigma = 1 $ . And it is declined with constant potential temperature (along dry adiabat line) below where air pressure is 1.0e+4. It is constant above where air pressure is 1.0e+4.
    • Specific humidity: 6.11641e-3 [kg kg-1] at $ sigma = 1 $ . And it is constant below where it is same as 75 % of saturation specific humidity. Above that, specific humidity is 75 % of saturation specific humidity.

Procedures List

SetInitData :初期値データ取得
———— :————
SetInitData :Get initial data

NAMELIST

NAMELIST#initial_data_nml

Methods

Included Modules

gridset dc_types dc_message axesset dc_string set_1d_profile sltt_debug constants saturate constants0 gauss_quad namelist_util dc_iounit

Public Instance methods

Subroutine :

This procedure input/output NAMELIST#initial_data_nml .

[Source]

  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

    ! 飽和比湿の算出
    ! Evaluate saturation specific humidity
    !
    use saturate, only: SaturateInit


    ! 宣言文 ; 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, Ueq, UGeos, VGeos
          !
          ! デフォルト値については初期化手続 "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


    ! デフォルト値の設定 (まずは 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
      QVapAvr = 0.0e0_DP
      Ueq     = 0.0_DP
      UGeos   = 0.0_DP
      VGeos   = 0.0_DP
    case ( 'agcm 5.3 default' )
      TempAvr = 250.0_DP
      PsAvr   = 1.0e+5_DP
      QVapAvr = 1.0e-10_DP
      Ueq     = 0.0_DP
      UGeos   = 0.0_DP
      VGeos   = 0.0_DP
    case ( 'sugiyama et al. (2008)' )
      TempAvr = 490.0_DP
      PsAvr   = 3.0e+6_DP
      QVapAvr = 6.11641e-3_DP
      Ueq     = 0.0_DP
      UGeos   = 0.0_DP
      VGeos   = 0.0_DP
    case ( 'polvani et al. (2004)' )
      TempAvr = 0.0_DP
      PsAvr   = 1.0e+5_DP
      QVapAvr = 0.0_DP
      Ueq     = 0.0_DP
      UGeos   = 0.0_DP
      VGeos   = 0.0_DP
    case ( 'venus' )
      TempAvr = 0.0_DP
      PsAvr   = 90.0e+5_DP
      QVapAvr = 0.0_DP
      Ueq     = 0.0_DP
      UGeos   = 0.0_DP
      VGeos   = 0.0_DP
    case ( '1-d profile' )
      TempAvr = 1.0e+100_DP
      PsAvr   = 1.0e+100_DP
      QVapAvr = 0.0_DP
      Ueq     = 0.0_DP
      UGeos   = 0.0_DP
      VGeos   = 0.0_DP
    case ( 'sltt debug' )
      TempAvr = 300.0_DP
      PsAvr   =   1.0e+5_DP
      QVapAvr =   0.0e0_DP
      Ueq     =  30.0_DP
      UGeos   =   0.0_DP
      VGeos   =   0.0_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


    ! Initialization of modules used in this module
    !
    ! 飽和比湿の算出
    ! Evaluate saturation specific humidity
    !
    call SaturateInit


    ! 印字 ; 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, '  Ueq     = %f', d = (/ Ueq  /) )
    call MessageNotify( 'M', module_name, '  UGeos   = %f', d = (/ UGeos /) )
    call MessageNotify( 'M', module_name, '  VGeos   = %f', d = (/ VGeos /) )
    call MessageNotify( 'M', module_name, '' )
    call MessageNotify( 'M', module_name, '-- version = %c', c1 = trim(version) )

    initial_data_inited = .true.

  end subroutine InitDataInit
Pattern
Variable :
Pattern :character(STRING), save, public
: 初期値データのパターン. 以下のパターンを選択可能.

Initial data pattern Available patterns are as follows.

  • "Small Disturbance of Temperature" (default value)
  • "AGCM 5.3 Default"
  • "Sugiyama et al. (2008)"
PsAvr
Variable :
PsAvr :real(DP), save, public
: $ bar{p_s} $ . 地表面気圧平均値. Mean surface pressure
QVapAvr
Variable :
QVapAvr :real(DP), save, public
: $ bar{q} $ . 比湿平均値. Mean specific humidity
Subroutine :
xyz_U(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(out)
: $ u $ . 東西風速. Eastward wind
xyz_V(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(out)
: $ v $ . 南北風速. Northward wind
xyz_Temp(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(out)
: $ T $ . 温度. Temperature
xyz_QVap(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(out)
: $ q $ . 比湿. Specific humidity
xy_Ps(0:imax-1, 1:jmax) :real(DP), intent(out)
: $ p_s $ . 地表面気圧. Surface pressure

初期値データのサンプルを提供します.

Prepare sample data of initial data

[Source]

  subroutine SetInitData( xyz_U, xyz_V, xyz_Temp, xyz_QVap, xy_Ps )
    !
    ! 初期値データのサンプルを提供します. 
    ! 
    ! Prepare sample data of initial data
    !

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

    ! 座標データ設定
    ! Axes data settings
    !
    use axesset, only: x_Lon, y_Lat, z_Sigma
                              ! $ \sigma $ レベル (整数). 
                              ! Full $ \sigma $ level

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

    ! ファイルから 1 次元プロファイルを読んで設定する. 
    ! read 1-D profile from a file and set it 
    !
    use set_1d_profile, only : Set1DProfilePs, Set1DProfileAtm

    ! セミラグ移流テスト用初期値作成
    ! Preparation of initial condition for SLTT advection
    !
    use sltt_debug, only : SLTTDebugSetUV, SLTTDebugSetQ


    ! 宣言文 ; Declaration statements
    !
    implicit none
    real(DP), intent(out):: xyz_U  (0:imax-1, 1:jmax, 1:kmax)
                              ! $ u $ .   東西風速. Eastward wind
    real(DP), intent(out):: xyz_V  (0:imax-1, 1:jmax, 1:kmax)
                              ! $ v $ .   南北風速. Northward wind
    real(DP), intent(out):: xyz_Temp  (0:imax-1, 1:jmax, 1:kmax)
                              ! $ T $ .   温度. Temperature
    real(DP), intent(out):: xyz_QVap  (0:imax-1, 1:jmax, 1:kmax)
                              ! $ q $ .   比湿. Specific humidity
    real(DP), intent(out):: xy_Ps (0:imax-1, 1:jmax)
                              ! $ p_s $ . 地表面気圧. Surface pressure

    ! 作業変数
    ! Work variables
    !
    real(DP) :: xyz_Press(0:imax-1, 1:jmax, 1:kmax)

    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 ) then
      call MessageNotify( 'E', module_name, 'This module has not been initialized.' )
    end if


    select case ( LChar( trim(Pattern) ) )
    case ( 'small disturbance of temperature' )
      !
      ! 微小な温度擾乱のある静止場
      ! Stationary field with small disturbance of temperature
      !

      xyz_U    = 0.0_DP
      xyz_V    = 0.0_DP
      xyz_Temp = TempAvr
      xy_Ps    = PsAvr
      xyz_QVap = QVapAvr

      ! 温度に擾乱を与える
      ! Add perturbation to temperature
      !
      do k = 1, kmax
        do j = 1, jmax
          do i = 0, imax - 1
            xyz_Temp(i,j,k) = xyz_Temp(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

      ! 東西風速を与える
      ! Add eastward wind
      !
      do j = 1, jmax
        xyz_U(:,j,:) = Ueq * cos(y_Lat(j))
      end do

    case ( 'agcm 5.3 default' )
      !
      ! AGCM5.3 のデフォルト初期値
      ! AGCM5.3 default initial values
      !

      xyz_U    = 0.0_DP
      xyz_V    = 0.0_DP
      xyz_Temp = TempAvr
      xy_Ps    = PsAvr
      xyz_QVap = QVapAvr

      ! 温度に擾乱を与える
      ! Add perturbation to temperature
      !
      do k = 1, kmax
        do j = 1, jmax
          do i = 0, imax - 1
            xyz_Temp(i,j,k) = xyz_Temp(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

      ! 東西風速を与える
      ! Add eastward wind
      !
      do j = 1, jmax
        xyz_U(:,j,:) = Ueq * cos(y_Lat(j))
      end do

    case ( 'sugiyama et al. (2008)' )
      !
      ! Sugiyama et al. (2008) の初期値
      ! Initial values of Sugiyama et al. (2008)
      !
      call Sugiyamaetal2008InitData( xyz_U, xyz_V, xyz_Temp, xyz_QVap, xy_Ps )

    case ( 'polvani et al. (2004)' )
      !
      ! Polvani et al. (2004) の初期値
      ! Initial values of Polvani et al. (2004)
      !
      call Polvanietal2004InitData( xyz_U, xyz_V, xyz_Temp, xyz_QVap, xy_Ps )

    case ( 'venus' )

      call VenusInitData( xyz_U, xyz_V, xyz_Temp, xyz_QVap, xy_Ps )

    case ( '1-d profile' )

      xyz_U = UGeos
      xyz_V = VGeos

      call Set1DProfilePs( xy_Ps )

      do k = 1, kmax
        xyz_Press(:,:,k) = xy_Ps * z_Sigma(k)
      end do
      call Set1DProfileAtm( xyz_Press, xyz_Temp, xyz_QVap )

    case ( 'sltt debug' )

      call SLTTDebugSetUV( Ueq, xyz_U, xyz_V )
      xyz_Temp = TempAvr
      xy_Ps    = PsAvr
      call SLTTDebugSetQ( xyz_QVap )

    end select

  end subroutine SetInitData
TempAvr
Variable :
TempAvr :real(DP), save, public
: $ bar{T} $ . 温度平均値. Mean temperature
UGeos
Variable :
UGeos :real(DP), save, public
: $ u_{g} $ . Eastward geostrophic wind
Ueq
Variable :
Ueq :real(DP), save, public
: $ u_{eq} $ . 赤道上の東西風速. Eastward wind on the equator
VGeos
Variable :
VGeos :real(DP), save, public
: $ v_{g} $ . Northward geostrophic wind
initial_data_inited
Variable :
initial_data_inited = .false. :logical, save, public
: 初期設定フラグ. Initialization flag

Private Instance methods

Subroutine :
xyz_U(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(out)
: $ u $ . 東西風速. Eastward wind
xyz_V(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(out)
: $ v $ . 南北風速. Northward wind
xyz_Temp(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(out)
: $ T $ . 温度. Temperature
xyz_QVap(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(out)
: $ q $ . 比湿. Specific humidity
xy_Ps(0:imax-1, 1:jmax) :real(DP), intent(out)
: $ p_s $ . 地表面気圧. Surface pressure

Polvani et al. (2004) の初期値

initial data by Polvani et al. (2004)

[Source]

  subroutine Polvanietal2004InitData( xyz_U, xyz_V, xyz_Temp, xyz_QVap, xy_Ps )
    !
    ! Polvani et al. (2004) の初期値
    !
    ! initial data by Polvani et al. (2004)
    !

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

    ! 物理・数学定数設定
    ! Physical and mathematical constants settings
    !
    use constants0, only: PI                    ! $ \pi $.
                              ! 円周率. Circular constant

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

    ! 座標データ設定
    ! Axes data settings
    !
    use axesset, only: x_Lon, y_Lat, z_Sigma, y_Lat_Weight
                              ! $ \Delta \varphi $ [rad.] .
                              ! 緯度座標重み.
                              ! Weight of latitude

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

    ! ガウス重み, 分点の計算
    ! Calculate Gauss node and Gaussian weight
    !
    use gauss_quad, only : GauLeg


    ! 宣言文 ; Declaration statements
    !
    implicit none
    real(DP), intent(out):: xyz_U  (0:imax-1, 1:jmax, 1:kmax)
                              ! $ u $ .   東西風速. Eastward wind
    real(DP), intent(out):: xyz_V  (0:imax-1, 1:jmax, 1:kmax)
                              ! $ v $ .   南北風速. Northward wind
    real(DP), intent(out):: xyz_Temp  (0:imax-1, 1:jmax, 1:kmax)
                              ! $ T $ .   温度. Temperature
    real(DP), intent(out):: xyz_QVap  (0:imax-1, 1:jmax, 1:kmax)
                              ! $ q $ .   比湿. Specific humidity
    real(DP), intent(out):: xy_Ps (0:imax-1, 1:jmax)
                              ! $ p_s $ . 地表面気圧. Surface pressure


    ! 作業変数
    ! Work variables
    !

    ! パラメタ設定
    ! 
    !
    !
    real(DP), parameter:: dz0 = 5.0d3
                              ! dz_0 [m].
    real(DP), parameter:: z0  = 22.0d3
                              ! z_0 [m]
    real(DP), parameter:: z1  = 30.0d3
                              ! z_1 [m]
    real(DP), parameter:: U0  = 50.0_DP
                              ! u_0 [m s^-1]

    real(DP), parameter:: ScaleHeight = 7.34d3
                              ! ScaleHeight [m]

    real(DP):: z_Height (1:kmax)
                              ! height [m]
    real(DP):: z_F (1:kmax)
                              ! function used to calculate zonal wind and temperature
    real(DP):: z_DFDZ (1:kmax)
                              ! z derivative of z_F
    real(DP):: z_TempUSStd(1:kmax)
                              ! Temperature by US Standard atmosphere [K]


    integer, parameter :: NGauQuad = 100
    real(DP)           :: a_GauQuadLat( 1:NGauQuad )
    real(DP)           :: a_GauWeight ( 1:NGauQuad )
    real(DP)           :: az_U        ( 1:NGauQuad, 1:kmax )
    real(DP)           :: az_DUDZ     ( 1:NGauQuad, 1:kmax )
    real(DP)           :: az_DTempDPhi( 1:NGauQuad, 1:kmax )
    real(DP)           :: z_Temp0     (1:kmax)


    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
    integer:: l               ! 緯度方向に回る DO ループ用作業変数
                              ! Work variables for DO loop in latitude


    ! 実行文 ; Executable statement

    if ( .not. initial_data_inited ) then
      call MessageNotify( 'E', module_name, 'This module has not been initialized.' )
    end if


    ! 高度の計算
    ! Calculate height
    !

    do k = 1, kmax
      z_Height = - ScaleHeight * log( z_Sigma )
    end do

    do k = 1, kmax
      if      ( z_Height(k) * 1.0d-3 < 11.0_DP ) then
        z_TempUSStd(k) = 288.15_DP - 6.5_DP * z_Height(k) * 1.0d-3
      else if ( z_Height(k) * 1.0d-3 < 20.0_DP ) then
        z_TempUSStd(k) = 216.65_DP
      else if ( z_Height(k) * 1.0d-3 < 32.0_DP ) then
        z_TempUSStd(k) = 216.65_DP + 1.0_DP * ( z_Height(k) * 1.0d-3 - 20.0_DP )
      else if ( z_Height(k) * 1.0d-3 < 47.0_DP ) then
        z_TempUSStd(k) = 228.65_DP + 2.8_DP * ( z_Height(k) * 1.0d-3 - 32.0_DP )
      else if ( z_Height(k) * 1.0d-3 < 51.0_DP ) then
        z_TempUSStd(k) = 270.65_DP
      else if ( z_Height(k) * 1.0d-3 < 71.0_DP ) then
        z_TempUSStd(k) = 270.65_DP - 2.8_DP * ( z_Height(k) * 1.0d-3 - 51.0_DP )
      else if ( z_Height(k) * 1.0d-3 < 80.0_DP ) then
        z_TempUSStd(k) = 214.65_DP - 2.0_DP * ( z_Height(k) * 1.0d-3 - 71.0_DP )
      else
        z_TempUSStd(k) = 196.65_DP
      end if
    end do

    z_F = 0.5d0 * ( 1.0d0 - tanh( ( z_Height - z0 ) / dz0 )**3 ) * sin( PI * z_Height / z1 )
    z_DFDZ = - 1.5_DP / dz0 * tanh( ( z_Height - z0 ) / dz0 )**2 / cosh( ( z_Height - z0 ) / dz0 )**2 * sin( PI * z_Height / z1 ) + 0.5_DP * ( 1.0_DP - tanh( ( z_Height - z0 ) / dz0 )**3 ) * PI / z1 * cos( PI * z_Height / z1 )



    xyz_U    = 0.0_DP
    xyz_V    = 0.0_DP
    xyz_Temp = 0.0_DP
    xy_Ps    = PsAvr
    xyz_QVap = QVapAvr


    ! 東西風速の計算
    ! Calculate eastward wind
    !
    do k = 1, kmax
      do j = 1, jmax
        if ( y_Lat(j) > 0.0_DP ) then
          xyz_U(:,j,k) = U0 * sin( PI * sin( y_Lat(j) )**2 )**3 * z_F(k)
        else
          xyz_U(:,j,k) =0.0_DP
        end if
      end do
    end do


    ! 温度の計算
    ! Calculate temperature
    !
    do j = 1, jmax

      call GauLeg( -PI/2.0_DP, y_Lat(j), NGauQuad, a_GauQuadLat, a_GauWeight )

      do k = 1, kmax
        do l = 1, NGauQuad
          if ( a_GauQuadLat(l) > 0.0_DP ) then
            az_U   (l,k) = U0 * sin( PI * sin( a_GauQuadLat(l) )**2 )**3 * z_F(k)
            az_DUDZ(l,k) = U0 * sin( PI * sin( a_GauQuadLat(l) )**2 )**3 * z_DFDZ(k)
          else
            az_U   (l,k) =0.0_DP
            az_DUDZ(l,k) =0.0_DP
          end if
        end do
      end do

      do k = 1, kmax
        do l = 1, NGauQuad
          if ( a_GauQuadLat(l) > 0.0_DP ) then
            az_DTempDPhi(l,k) = - ScaleHeight / GasRDry * (  2.0_DP * RPlanet * Omega * sin( a_GauQuadLat(l) ) + 2.0_DP * az_U(l,k) * tan( a_GauQuadLat(l) ) ) * az_DUDZ(l,k)
          else
            az_DTempDPhi(l,k) = 0.0_DP
          end if
        end do
      end do

      do k = 1, kmax
        xyz_Temp(:,j,k) = 0.0_DP
      end do
      do k = 1, kmax
        do l = 1, NGauQuad
          xyz_Temp(:,j,k) = xyz_Temp(:,j,k) + az_DTempDPhi(l,k) * a_GauWeight(l)
        end do
      end do

    end do

    !  Calculate T0
    !
    do k = 1, kmax
      z_Temp0(k) = 0.0_DP
      do j = 1, jmax
        z_Temp0(k) = z_Temp0(k) + xyz_Temp(0,j,k) * y_Lat_weight(j)
      end do
    end do
    z_Temp0 = z_Temp0 / 2.0_DP
    !
    !  add T0
    !
    do k = 1, kmax
      xyz_Temp(:,:,k) = xyz_Temp(:,:,k) - z_Temp0(k) + z_TempUSStd(k)
    end do


    ! Code for debug
!!$    do k = 1, kmax
!!$      z_Temp0(k) = 0.0_DP
!!$      do j = 1, jmax
!!$        z_Temp0(k) = z_Temp0(k) + xyz_Temp(0,j,k) * y_Lat_weight(j)
!!$      end do
!!$    end do
!!$    z_Temp0 = z_Temp0 / 2.0_DP
!!$    do k = 1, kmax
!!$      write( 6, * ) k, z_Temp0(k), z_TempUSStd(k), z_Temp0(k)-z_TempUSStd(k)
!!$    end do
!!$    stop



    ! 温度に擾乱を与える
    ! Add perturbation to temperature
    !
    do k = 1, kmax
      do j = 1, jmax
        do i = 0, imax-1
          if ( ( x_Lon(i) >= 0 ) .and. ( x_Lon(i) < PI ) ) then
            xyz_Temp(i,j,k) = xyz_Temp(i,j,k) + 1.0_DP / cosh( 3.0_DP * ( x_Lon(i)               ) )**2 * 1.0_DP / cosh( 6.0_DP * ( y_Lat(j) - PI / 4.0_DP ) )**2
          else if ( ( x_Lon(i) >= PI ) .and. ( x_Lon(i) < 2.0_DP * PI ) ) then
            xyz_Temp(i,j,k) = xyz_Temp(i,j,k) + 1.0_DP / cosh( 3.0_DP * ( x_Lon(i) - 2.0_DP * PI ) )**2 * 1.0_DP / cosh( 6.0_DP * ( y_Lat(j) - PI / 4.0_DP ) )**2
          end if
        end do
      end do
    end do


  end subroutine Polvanietal2004InitData
Subroutine :
xyz_U(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(out)
: $ u $ . 東西風速. Eastward wind
xyz_V(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(out)
: $ v $ . 南北風速. Northward wind
xyz_Temp(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(out)
: $ T $ . 温度. Temperature
xyz_QVap(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(out)
: $ q $ . 比湿. Specific humidity
xy_Ps(0:imax-1, 1:jmax) :real(DP), intent(out)
: $ p_s $ . 地表面気圧. Surface pressure

Sugiyama et al. (2008) の初期値 Initial values of Sugiyama et al. (2008)

[Source]

  subroutine Sugiyamaetal2008InitData( xyz_U, xyz_V, xyz_Temp, xyz_QVap, xy_Ps )
    !
    ! Sugiyama et al. (2008) の初期値
    ! Initial values of Sugiyama et al. (2008)
    !

    ! モジュール引用 ; 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

    ! 飽和比湿の算出
    ! Evaluate saturation specific humidity
    !
    use saturate, only: xyz_CalcQVapSat

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

    ! 宣言文 ; Declaration statements
    !
    implicit none
    real(DP), intent(out):: xyz_U  (0:imax-1, 1:jmax, 1:kmax)
                              ! $ u $ .   東西風速. Eastward wind
    real(DP), intent(out):: xyz_V  (0:imax-1, 1:jmax, 1:kmax)
                              ! $ v $ .   南北風速. Northward wind
    real(DP), intent(out):: xyz_Temp  (0:imax-1, 1:jmax, 1:kmax)
                              ! $ T $ .   温度. Temperature
    real(DP), intent(out):: xyz_QVap  (0:imax-1, 1:jmax, 1:kmax)
                              ! $ q $ .   比湿. Specific humidity
    real(DP), intent(out):: xy_Ps (0:imax-1, 1:jmax)
                              ! $ p_s $ . 地表面気圧. 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 ) then
      call MessageNotify( 'E', module_name, 'This module has not been initialized.' )
    end if


    xyz_U    = 0.0_DP
    xyz_V    = 0.0_DP
    xyz_Temp = TempAvr
    xy_Ps    = PsAvr
    xyz_QVap = QVapAvr

    ! 温度の計算
    ! Calculate temperature
    !
    xyz_PotTemp = TempAvr
    xy_TempMin  = TempAvr

    do k = 1, kmax
      xyz_Temp(:,:,k) = xyz_PotTemp(:,:,k) * ( z_Sigma(k) )**( GasRDry / CpDry )

      if ( PsAvr * z_Sigma(k) < 1.0e+4_DP ) then
        xyz_Temp(:,:,k) = xy_TempMin
      else
        xy_TempMin = xyz_Temp(:,:,k)
      end if
    end do

    ! 温度に擾乱を与える
    ! Add perturbation to temperature
    !
    do k = 1, kmax
      do j = 1, jmax
        do i = 0, imax - 1
          xyz_Temp(i,j,k) = xyz_Temp(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_Ps * z_Sigma(k)
    end do

    xyz_QVapSat = xyz_CalcQVapSat( xyz_Temp, xyz_Press )

    ! 比湿の計算
    ! Calculate specific humidity
    !
    where ( xyz_QVap > xyz_QVapSat * 0.75 )
      xyz_QVap = xyz_QVapSat * 0.75
    end where


  end subroutine Sugiyamaetal2008InitData
Subroutine :
xyz_U(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(out)
: $ u $ . 東西風速. Eastward wind
xyz_V(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(out)
: $ v $ . 南北風速. Northward wind
xyz_Temp(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(out)
: $ T $ . 温度. Temperature
xyz_QVap(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(out)
: $ q $ . 比湿. Specific humidity
xy_Ps(0:imax-1, 1:jmax) :real(DP), intent(out)
: $ p_s $ . 地表面気圧. Surface pressure

[Source]

  subroutine VenusInitData( xyz_U, xyz_V, xyz_Temp, xyz_QVap, xy_Ps )

    ! 種別型パラメタ
    ! Kind type parameter
    !
    use dc_types, only: DP, STRING     ! 文字列.       Strings.


    ! 格子点設定
    ! Grid points settings
    !
    use gridset, only: imax, jmax, kmax    ! 鉛直層数.
                               ! Number of vertical level

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

    use axesset, only: y_Lat, z_Sigma, r_Sigma, r_DelSigma            ! $ \Delta \sigma $ (半整数).
                              ! $ \Delta \sigma $ (half)

    real(DP), intent(out):: xyz_U  (0:imax-1, 1:jmax, 1:kmax)
                              ! $ u $ .   東西風速. Eastward wind
    real(DP), intent(out):: xyz_V  (0:imax-1, 1:jmax, 1:kmax)
                              ! $ v $ .   南北風速. Northward wind
    real(DP), intent(out):: xyz_Temp  (0:imax-1, 1:jmax, 1:kmax)
                              ! $ T $ .   温度. Temperature
    real(DP), intent(out):: xyz_QVap  (0:imax-1, 1:jmax, 1:kmax)
                              ! $ q $ .   比湿. Specific humidity
    real(DP), intent(out):: xy_Ps (0:imax-1, 1:jmax)
                              ! $ p_s $ . 地表面気圧. Surface pressure


    !
    ! local variables
    !
    real(DP) :: SurfTemp
    real(DP) :: xyr_Temp     (0:imax-1,1:jmax,0:kmax)
    real(DP) :: xy_SurfHeight(0:imax-1,1:jmax)
    real(DP) :: xyz_Height   (0:imax-1,1:jmax,1:kmax)
    real(DP) :: z( 5 ), a( 6 ), ah( 5 ), d( 5 )
    integer  :: j
    integer  :: k
    integer  :: l
    integer  :: m


    ! Coefficients for thermal structure by Hou and Farrel (1987)
    !
    z ( 1 ) =   0.0d3
    z ( 2 ) =  10.0d3
    z ( 3 ) =  25.0d3
    z ( 4 ) =  55.0d3
    z ( 5 ) = 100.0d3

    ah( 1 ) =  -1.0d-3
    ah( 2 ) =  -1.0d-3
    ah( 3 ) =  -3.1d-3
    ah( 4 ) =  -6.75d-3
    ah( 5 ) =  10.0d-3

    d ( 1 ) =  10.0d3
    d ( 2 ) =  10.0d3
    d ( 3 ) =   8.0d3
    d ( 4 ) =   5.0d3
    d ( 5 ) =  70.0d3


    a ( 1 ) =   0.0d0
    do l = 2, 6
      a( l ) = 2.0d0 * ah( l-1 ) * d( l-1 ) + a( l-1 )
    end do


    SurfTemp      = 750.0_DP
    xy_SurfHeight =   0.0_DP


    ! Initialization
    xyz_Temp = 200.0_DP

    ! Iteration
    do m = 1, 10

      xyr_Temp(:,:,0) = xyz_Temp(:,:,1)
      do k = 1, kmax-1
        xyr_Temp(:,:,k) = ( xyz_Temp(:,:,k) + xyz_Temp(:,:,k+1) ) / 2.0_DP
      end do
      xyr_Temp(:,:,kmax) = xyz_Temp(:,:,kmax)


      xyz_Height(:,:,1) = xy_SurfHeight + GasRDry / Grav * xyz_Temp(:,:,1) * ( 1. - z_Sigma(1) )
      do k = 2, kmax
        xyz_Height(:,:,k) = xyz_Height(:,:,k-1) + GasRDry / Grav * xyr_Temp(:,:,k-1) * r_DelSigma(k-1) / r_Sigma(k-1)
      end do

      xyz_Temp = SurfTemp - Grav / CpDry * xyz_Height
      do l = 1, 5
        xyz_Temp = xyz_Temp - ( a(l+1) - a(l) ) * 0.5d0 * ( 1.0d0 + tanh( ( 0.0d0      - z(l) ) / d(l) ) )
        xyz_Temp = xyz_Temp + ( a(l+1) - a(l) ) * 0.5d0 * ( 1.0d0 + tanh( ( xyz_Height - z(l) ) / d(l) ) )
      end do

    end do

    ! add perturbation
    xyz_Temp(0,1,1) = xyz_Temp(0,1,1) + 1.0_DP


    do k = 1, kmax
      do j = 1, jmax
        xyz_U(:,j,k) = Ueq * cos(y_Lat(j))
      end do
    end do
    xyz_V    = 0.0_DP
    xyz_QVap = QVapAvr
    xy_Ps    = PsAvr


  end subroutine VenusInitData
module_name
Constant :
module_name = ‘initial_data :character(*), parameter
: モジュールの名称. Module name
version
Constant :
version = ’$Name: dcpam5-20140314 $’ // ’$Id: initial_data.f90,v 1.13 2014-02-14 08:12:21 yot Exp $’ :character(*), parameter
: モジュールのバージョン Module version