Class phy_sponge_layer
In: physics/phy_sponge_layer.f90

速度を減衰させるためのスポンジ層

Sponge layer for damping velocity

Note that Japanese and English are described in parallel.

ある一定気圧 (デフォルトでは 1.0e4 Pa) 以下の層で, 速度の擾乱成分を潰すための減衰の効果を算出します. 減衰の強さは, ある時定数 (デフォルトは 1 日) で東西平均値に合わせるように決めます.

Calculate effect of damping perturbation of velocity to layers under certain pressure (by default, 1.0e4 Pa). Effect of damping is decided to Strength of damping is decided to match to east and west mean value by a certain time constant (default value is 1 day).

Procedures List

PhySpoLayCreate :PHYSPOLAY 型変数の初期設定
Damping :減衰の効果の算出
PhySpoLayClose :PHYSPOLAY 型変数の終了処理
PhySpoLayPutLine :PHYSPOLAY 型変数に格納されている情報の印字
PhySpoLayinitialized :PHYSPOLAY 型変数が初期設定されているか否か
PhySpoLaySetTime :時刻の設定
———— :————
PhySpoLayCreate :Constructor of "PHYSPOLAY"
Damping :Calculation effect of damping
PhySpoLayClose :Deconstructor of "PHYSPOLAY"
PhySpoLayPutLine :Print information of "PHYSPOLAY"
PhySpoLayinitialized :Check initialization of "PHYSPOLAY"
PhySpoLaySetTime :Configure time

Usage

始めに, PHYSPOLAY 型の変数を定義し, PhySpoLayCreate で初期設定を行います. Damping に東西風速, 南北風速, 地表面気圧を与えることで, それぞれへの減衰の効果を算出します. PHYSPOLAY 型の変数の終了処理には PhySpoLayClose を用いてください.

First, initialize "PHYSPOLAY" by "PhySpoLayCreate". "Damping" receives zonal and meridional wind and surface pressure, and returns effect of dumping them. In order to terminate "PHYSPOLAY", use "PhySpoLayClose".

Methods

Included Modules

dc_types dc_date_types gt4_history_nmlinfo intavr_operate dc_trace dc_string dc_present dc_message dcpam_error dc_error dc_date dc_hash gt4_history dc_iounit

Public Instance methods

Damping( phy_spo_lay, xyz_U, xyz_V, xy_Ps, xyz_DUDtSpoDamping, xyz_DVDtSpoDamping, [historyput_flag], [err] )
Subroutine :
phy_spo_lay :type(PHYSPOLAY), intent(inout)
xyz_U(0:phy_spo_lay%imax-1, 0:phy_spo_lay%jmax-1, 0:phy_spo_lay%kmax-1) :real(DP), intent(in)
: $ u $ . 東西風速. Zonal wind
xyz_V(0:phy_spo_lay%imax-1, 0:phy_spo_lay%jmax-1, 0:phy_spo_lay%kmax-1) :real(DP), intent(in)
: $ v $ . 南北風速. Meridional wind
xy_Ps(0:phy_spo_lay%imax-1, 0:phy_spo_lay%jmax-1) :real(DP), intent(in)
: $ p_s $ . 地表面気圧. Surface pressure
xyz_DUDtSpoDamping(0:phy_spo_lay%imax-1, 0:phy_spo_lay%jmax-1, 0:phy_spo_lay%kmax-1) :real(DP), intent(out)
: $ DP{u}{t} $ . 減衰の効果による東西風速変化. Zonal wind tendency by damping effect
xyz_DVDtSpoDamping(0:phy_spo_lay%imax-1, 0:phy_spo_lay%jmax-1, 0:phy_spo_lay%kmax-1) :real(DP), intent(out)
: $ DP{v}{t} $ . 減衰の効果による南北風速変化. Meridional wind tendency by damping effect
historyput_flag :logical, intent(in), optional
: データ出力のフラグ. SetTime によって時刻を明示的に 指定した場合には, この引数に .true. または .false. を指定する ことでデータ出力のオンオフを 明示的に指定する必要があります. デフォルトは .false. です.

Data output flag. When time is specified by "SetTime", explicit specification of data output on/off by specifying ".true." or ".false." to this argument. Default value is ".false.".

err :logical, intent(out), optional
: 例外処理用フラグ. デフォルトでは, この手続き内でエラーが 生じた場合, プログラムは強制終了します. 引数 err が与えられる場合, プログラムは強制終了せず, 代わりに err に .true. が代入されます.

Exception handling flag. By default, when error occur in this procedure, the program aborts. If this err argument is given, .true. is substituted to err and the program does not abort.

東西風速 xyz_U, 南北風速 xyz_V および 地表面気圧 xy_Ps から, 減衰の効果による東西風速と南北風速の変化率 xyz_DUDtSpoDamping, xyz_DVDtSpoDamping を算出して返します. 計算後, phy_spo_lay 内の時刻を更新します.

なお, 与えられた phy_spo_layPhySpoLayCreate によって初期設定されていない場合, プログラムはエラーを発生させます.

Calculate and return zonal wind and meridional tendency by damping in sponge layer xyz_DUDtSpoDamping, xyz_DVDtSpoDamping from zonal wind xyz_U and meridional wind xyz_V and surface pressure xy_Ps. After the calculation, time stored in phy_spo_lay is updated.

If phy_spo_lay is not initialized by "PhySpoLayCreate" yet, error is occurred.

Alias for PhySpoLayDamping

PHYSPOLAY
Derived Type :
initialized = .false. :logical
: 初期設定フラグ. Initialization flag
imax :integer
: 経度格子点数. Number of grid points in longitude
jmax :integer
: 緯度格子点数. Number of grid points in latitude
kmax :integer
: 鉛直層数. Number of vertical level
x_Lon(:) =>null() :real(DP), pointer
: 経度. Longitude
x_Lon_Weight(:) =>null() :real(DP), pointer
: 経度積分用座標重み. Weight for integration in longitude
y_Lat(:) =>null() :real(DP), pointer
: 緯度. Latitude
y_Lat_Weight(:) =>null() :real(DP), pointer
: 緯度積分用座標重み. Weight for integration in latitude
z_Sigma(:) =>null() :real(DP), pointer
: $ sigma $ レベル (整数). Full $ sigma $ level
z_DelSigma(:) =>null() :real(DP), pointer
: $ Delta sigma $ (整数). $ Delta sigma $ (Full)
DampingTime = 86400.0_DP :real(DP)
: 減衰の時定数 (単位: 秒). Damping time constant (units: seconds)
UpperPress = 0.0_DP :real(DP)
: 減衰範囲の上端の気圧. Pressure in the top within the range of damping
LowerPress = 1.0e4_DP :real(DP)
: 減衰範囲の下端の気圧. Pressure at bottom within the range of damping
current_time :type(DC_DIFFTIME)
: 現在時刻. Current time.
delta_time :type(DC_DIFFTIME)
: $ Delta t $ . タイムステップ. Time step
intavr_opr :type(INTAVROPR)
: intavr_operate モジュールのオブジェクト. Object of "intavr_operate" module
gthstnml :type(GTHST_NMLINFO)
: NAMELIST#phy_sponge_layer_history_nml から入手される個別のデータ出力情報.

Individual data output information from "NAMELIST#phy_sponge_layer_history_nml".

まず, PhySpoLayCreate で "PHYSPOLAY" 型の 変数を初期設定して下さい. 初期設定された "PHYSPOLAY" 型の変数を再度利用する際には, PhySpoLayClose によって終了処理を行ってください.

Initialize "PHYSPOLAY" variable by "PhySpoLayCreate" before usage. If you reuse "PHYSPOLAY" variable again for another application, terminate by "PhySpoLayClose".

Subroutine :
phy_spo_lay :type(PHYSPOLAY), intent(inout)
err :logical, intent(out), optional
: 例外処理用フラグ. デフォルトでは, この手続き内でエラーが 生じた場合, プログラムは強制終了します. 引数 err が与えられる場合, プログラムは強制終了せず, 代わりに err に .true. が代入されます.

Exception handling flag. By default, when error occur in this procedure, the program aborts. If this err argument is given, .true. is substituted to err and the program does not abort.

PHYSPOLAY 型の変数の終了処理を行います. なお, 与えられた phy_spo_layPhySpoLayCreate によって初期設定されていない場合, プログラムはエラーを発生させます.

Deconstructor of "PHYSPOLAY". Note that if phy_spo_lay is not initialized by "PhySpoLayCreate" yet, error is occurred.

[Source]

  subroutine PhySpoLayClose( phy_spo_lay, err )
    !
    ! PHYSPOLAY 型の変数の終了処理を行います. 
    ! なお, 与えられた *phy_spo_lay* が PhySpoLayCreate 
    ! によって初期設定されていない場合, プログラムはエラーを発生させます. 
    !
    ! Deconstructor of "PHYSPOLAY". 
    ! Note that if *phy_spo_lay* is not initialized 
    ! by "PhySpoLayCreate" yet, error is occurred. 
    !
    use dc_trace, only: BeginSub, EndSub
    use dc_string, only: PutLine, Printf, Split, StrInclude, StoA, JoinChar
    use dc_types, only: DP, STRING, TOKEN, STDOUT
    use dc_error, only: StoreError, DC_NOERR, DC_ENOTINIT
    use gt4_history_nmlinfo, only: HstNmlInfoClose, HstNmlInfoNames, HstNmlInfoAssocGtHist, HstNmlInfoPutLine
    use gt4_history, only: GT_HISTORY, HistoryClose, HistoryInitialized
    use intavr_operate, only: IntAvrOprClose
    implicit none
    type(PHYSPOLAY), intent(inout):: phy_spo_lay
    logical, intent(out), optional:: err
                              ! 例外処理用フラグ. 
                              ! デフォルトでは, この手続き内でエラーが
                              ! 生じた場合, プログラムは強制終了します. 
                              ! 引数 *err* が与えられる場合, 
                              ! プログラムは強制終了せず, 代わりに
                              ! *err* に .true. が代入されます. 
                              !
                              ! Exception handling flag. 
                              ! By default, when error occur in 
                              ! this procedure, the program aborts. 
                              ! If this *err* argument is given, 
                              ! .true. is substituted to *err* and 
                              ! the program does not abort. 

    !-----------------------------------
    !  ヒストリファイルへのデータ出力設定
    !  Configure the settings for history data output
    character(STRING):: name = ''
                              ! 変数名. Variable identifier
    character(STRING):: varnames
                              ! 変数名リスト. 
                              ! List of variables
    character(TOKEN), pointer:: varnames_array(:) =>null()
                              ! 変数名リスト配列. 
                              ! List of variables (array) 
    integer:: i, vnmax
    type(GT_HISTORY), pointer:: gthist =>null()
                              ! gt4_history モジュール用構造体. 
                              ! Derived type for "gt4_history" module

    !-----------------------------------
    !  作業変数
    !  Work variables
    integer:: stat
    character(STRING):: cause_c
    character(*), parameter:: subname = 'PhySpoLayClose'
  continue
    call BeginSub( subname )
    stat = DC_NOERR
    cause_c = ''

    !-----------------------------------------------------------------
    !  初期設定のチェック
    !  Check initialization
    !-----------------------------------------------------------------
    if ( .not. phy_spo_lay % initialized ) then
      stat = DC_ENOTINIT
      cause_c = 'PHYSPOLAY'
      goto 999
    end if

    !-----------------------------------------------------------------
    !  "PHYSPOLAY" の設定の消去
    !  Clear the settings for "PHYSPOLAY"
    !-----------------------------------------------------------------
    deallocate( phy_spo_lay % x_Lon )
    deallocate( phy_spo_lay % x_Lon_Weight )
    deallocate( phy_spo_lay % y_Lat )
    deallocate( phy_spo_lay % y_Lat_Weight )
    deallocate( phy_spo_lay % z_Sigma )
    deallocate( phy_spo_lay % z_DelSigma )
    phy_spo_lay % DampingTime = 86400.0_DP
    phy_spo_lay % UpperPress = 0.0_DP
    phy_spo_lay % LowerPress = 1.0e4_DP

    call IntAvrOprClose( intavr_opr = phy_spo_lay % intavr_opr )

    !-----------------------------------------------------------------
    !  ヒストリファイルへのデータ出力の終了処理
    !  Terminate the settings for history data output
    !-----------------------------------------------------------------
    varnames = HstNmlInfoNames( phy_spo_lay % gthstnml )
    call Split( str = varnames, sep = ',', carray = varnames_array )            ! (out)
    vnmax = size( varnames_array )

    do i = 1, vnmax
      name = varnames_array(i)
      if ( trim( name ) == '' ) exit
      nullify( gthist )
      call HstNmlInfoAssocGtHist( gthstnml = phy_spo_lay % gthstnml, name = name, history = gthist, err = err )                              ! (out)
      if ( HistoryInitialized( gthist ) ) then

        call HistoryClose( history = gthist, err = err )                        ! (out)

      end if
    end do

    !-----------------------------------------------------------------
    !  ヒストリファイルへのデータ出力設定の割付解除
    !  Deallocate the settings for history data output
    !-----------------------------------------------------------------
    call HstNmlInfoClose( phy_spo_lay % gthstnml, err = err )                   ! (out)

    !-----------------------------------------------------------------
    !  終了処理, 例外処理
    !  Termination and Exception handling
    !-----------------------------------------------------------------
    phy_spo_lay % initialized = .false.
999 continue
    call StoreError( stat, subname, err, cause_c )
    call EndSub( subname )

  end subroutine PhySpoLayClose
Subroutine :
phy_spo_lay :type(PHYSPOLAY), intent(inout)
imax :integer, intent(in)
: 経度格子点数. Number of grid points in longitude
jmax :integer, intent(in)
: 緯度格子点数. Number of grid points in latitude
kmax :integer, intent(in)
: 鉛直層数. Number of vertical level
x_Lon(0:imax-1) :real(DP), intent(in)
: 経度. Longitude
y_Lat(0:jmax-1) :real(DP), intent(in)
: 緯度. Latitude
z_Sigma(0:kmax-1) :real(DP), intent(in)
: $ sigma $ レベル (整数). Full $ sigma $ level
DelTime :real(DP), intent(in)
: $ Delta t $ . タイムステップ. Time step
DampingTime :real(DP), intent(in), optional
: 減衰の時定数 (単位: 秒). Damping time constant (units: seconds)
UpperPress :real(DP), intent(in), optional
: 減衰範囲の上端の気圧. Pressure in the top within the range of damping
LowerPress :real(DP), intent(in), optional
: 減衰範囲の下端の気圧. Pressure at bottom within the range of damping
x_Lon_Weight(0:imax-1) :real(DP), intent(in), optional
: 経度積分用座標重み. Weight for integration in longitude
y_Lat_Weight(0:jmax-1) :real(DP), intent(in), optional
: 緯度積分用座標重み. Weight for integration in latitude
z_DelSigma(0:kmax-1) :real(DP), intent(in), optional
: $ Delta sigma $ (整数). $ Delta sigma $ (Full)
current_time_value :real, intent(in), optional
: 現在時刻の数値. Numerical value of current time
current_time_unit :character(*), intent(in), optional
: 現在時刻の単位. Unit of current time
history_varlist :character(*), intent(in), optional
: ヒストリデータの出力変数リスト. カンマで区切って並べる. (例: "Data1,Data2" ).

List of variables output to history data. Delimiter is comma. (exp. "Data1,Data2" ).

history_interval_value :real, intent(in), optional
: ヒストリデータの出力間隔の数値. Numerical value for interval of history data output
history_interval_unit :character(*), intent(in), optional
: ヒストリデータの出力間隔の単位. Unit for interval of history data output
history_precision :character(*), intent(in), optional
: ヒストリデータの精度. Precision of history data
history_fileprefix :character(*), intent(in), optional
: ヒストリデータのファイル名の接頭詞. Prefix of history data filenames
nmlfile :character(*), intent(in), optional
: NAMELIST ファイルの名称. この引数に空文字以外を与えた場合, 指定されたファイルから NAMELIST 変数群を読み込みます. ファイルを読み込めない場合にはエラーを 生じます.

NAMELIST 変数群の詳細に関しては NAMELIST#phy_sponge_layer_nml を参照してください.

NAMELIST file name. If nonnull character is specified to this argument, NAMELIST group name is loaded from the file. If the file can not be read, an error occurs.

See "NAMELIST#phy_sponge_layer_nml" for details about a NAMELIST group.

err :logical, intent(out), optional
: 例外処理用フラグ. デフォルトでは, この手続き内でエラーが 生じた場合, プログラムは強制終了します. 引数 err が与えられる場合, プログラムは強制終了せず, 代わりに err に .true. が代入されます.

Exception handling flag. By default, when error occur in this procedure, the program aborts. If this err argument is given, .true. is substituted to err and the program does not abort.

PHYSPOLAY 型の変数の初期設定を行います. 他のサブルーチンを使用する前に必ずこのサブルーチンによって PHYSPOLAY 型の変数を初期設定してください.

なお, 与えられた phy_spo_lay が既に初期設定されている場合, プログラムはエラーを発生させます.

NAMELIST を利用する場合には引数 nmlfileNAMELIST ファイル名 を与えてください. NAMELIST 変数群の詳細に関しては NAMELIST#phy_sponge_layer_nml を参照してください.

Constructor of "PHYSPOLAY". Initialize phy_spo_lay by this subroutine, before other procedures are used,

Note that if phy_spo_lay is already initialized by this procedure, error is occurred.

In order to use NAMELIST, specify a NAMELIST filename to argument nmlfile. See "NAMELIST#phy_sponge_layer_nml" for details about a NAMELIST group.

[Source]

  subroutine PhySpoLayCreate( phy_spo_lay, imax, jmax, kmax, x_Lon, y_Lat, z_Sigma, DelTime, DampingTime, UpperPress, LowerPress, x_Lon_Weight, y_Lat_Weight, z_DelSigma, current_time_value, current_time_unit, history_varlist, history_interval_value, history_interval_unit, history_precision, history_fileprefix, nmlfile, err )
    !
    ! PHYSPOLAY 型の変数の初期設定を行います. 
    ! 他のサブルーチンを使用する前に必ずこのサブルーチンによって 
    ! PHYSPOLAY 型の変数を初期設定してください. 
    !
    ! なお, 与えられた *phy_spo_lay* が既に初期設定されている場合, 
    ! プログラムはエラーを発生させます. 
    !
    ! NAMELIST を利用する場合には引数 *nmlfile* に NAMELIST ファイル名
    ! を与えてください. NAMELIST 変数群の詳細に関しては 
    ! NAMELIST#phy_sponge_layer_nml を参照してください. 
    !
    ! Constructor of "PHYSPOLAY". 
    ! Initialize *phy_spo_lay* by this subroutine, 
    ! before other procedures are used, 
    !
    ! Note that if *phy_spo_lay* is already initialized 
    ! by this procedure, error is occurred. 
    !
    ! In order to use NAMELIST, specify a NAMELIST filename to 
    ! argument *nmlfile*. See "NAMELIST#phy_sponge_layer_nml" 
    ! for details about a NAMELIST group. 
    !
    use dc_trace, only: BeginSub, EndSub
    use dc_string, only: PutLine, Printf, Split, StrInclude, StoA, JoinChar
    use dc_types, only: DP, STRING, TOKEN, STDOUT
    use dc_present, only: present_and_not_empty, present_and_true, present_select
    use dc_message, only: MessageNotify
    use dcpam_error, only: StoreError, DCPAM_ESMALLVAL
    use dc_error, only: DC_NOERR, DC_EALREADYINIT, DC_EARGLACK, DC_ENEGATIVE, DC_ENOFILEREAD, USR_ERRNO, HST_EBADVARNAME
    use dc_date, only: DCDiffTimeCreate
    use dc_hash, only: HASH, DCHashPut, DCHashRewind, DCHashNext, DCHashDelete
    use gt4_history_nmlinfo, only: HstNmlInfoCreate, HstNmlInfoAdd, HstNmlInfoEndDefine, HstNmlInfoPutLine, HstNmlInfoAllNameValid
    use gt4_history, only: GT_HISTORY, HistoryAddVariable, HistoryAddAttr
    use intavr_operate, only: IntAvrOprCreate
    implicit none
    type(PHYSPOLAY), intent(inout):: phy_spo_lay
    integer, intent(in):: imax ! 経度格子点数. 
                               ! Number of grid points in longitude
    integer, intent(in):: jmax ! 緯度格子点数. 
                               ! Number of grid points in latitude
    integer, intent(in):: kmax ! 鉛直層数. 
                               ! Number of vertical level
    real(DP), intent(in):: x_Lon (0:imax-1)
                              ! 経度. Longitude
    real(DP), intent(in):: y_Lat (0:jmax-1)
                              ! 緯度. Latitude
    real(DP), intent(in):: z_Sigma (0:kmax-1)
                              ! $ \sigma $ レベル (整数). 
                              ! Full $ \sigma $ level
    real(DP), intent(in):: DelTime    ! $ \Delta t $ . タイムステップ. Time step
    real(DP), intent(in), optional:: DampingTime
                              ! 減衰の時定数 (単位: 秒). 
                              ! Damping time constant (units: seconds)
    real(DP), intent(in), optional:: UpperPress
                              ! 減衰範囲の上端の気圧. 
                              ! Pressure in the top within the range of damping
    real(DP), intent(in), optional:: LowerPress
                              ! 減衰範囲の下端の気圧. 
                              ! Pressure at bottom within the range of damping
    real(DP), intent(in), optional:: x_Lon_Weight (0:imax-1)
                              ! 経度積分用座標重み. 
                              ! Weight for integration in longitude
    real(DP), intent(in), optional:: y_Lat_Weight (0:jmax-1)
                              ! 緯度積分用座標重み. 
                              ! Weight for integration in latitude
    real(DP), intent(in), optional:: z_DelSigma (0:kmax-1)
                              ! $ \Delta \sigma $ (整数). 
                              ! $ \Delta \sigma $ (Full)
    real, intent(in), optional:: current_time_value
                              ! 現在時刻の数値. Numerical value of current time
    character(*), intent(in), optional:: current_time_unit
                              ! 現在時刻の単位. Unit of current time
    character(*), intent(in), optional:: history_varlist
                              ! ヒストリデータの出力変数リスト. 
                              ! カンマで区切って並べる. 
                              ! (例: "Data1,Data2" ). 
                              ! 
                              ! List of variables output to history data. 
                              ! Delimiter is comma. 
                              ! (exp. "Data1,Data2" ). 
                              ! 
    real, intent(in), optional:: history_interval_value
                              ! ヒストリデータの出力間隔の数値. 
                              ! Numerical value for interval of history data output
    character(*), intent(in), optional:: history_interval_unit
                              ! ヒストリデータの出力間隔の単位. 
                              ! Unit for interval of history data output
    character(*), intent(in), optional:: history_precision
                              ! ヒストリデータの精度. 
                              ! Precision of history data
    character(*), intent(in), optional:: history_fileprefix
                              ! ヒストリデータのファイル名の接頭詞. 
                              ! Prefix of history data filenames

    character(*), intent(in), optional:: nmlfile
                              ! NAMELIST ファイルの名称. 
                              ! この引数に空文字以外を与えた場合, 
                              ! 指定されたファイルから 
                              ! NAMELIST 変数群を読み込みます. 
                              ! ファイルを読み込めない場合にはエラーを
                              ! 生じます. 
                              !
                              ! NAMELIST 変数群の詳細に関しては 
                              ! NAMELIST#phy_sponge_layer_nml 
                              ! を参照してください. 
                              !
                              ! NAMELIST file name. 
                              ! If nonnull character is specified to 
                              ! this argument, 
                              ! NAMELIST group name is loaded from the 
                              ! file. 
                              ! If the file can not be read, 
                              ! an error occurs.
                              ! 
                              ! See "NAMELIST#phy_sponge_layer_nml" 
                              ! for details about a NAMELIST group. 
                              ! 

    logical, intent(out), optional:: err
                              ! 例外処理用フラグ. 
                              ! デフォルトでは, この手続き内でエラーが
                              ! 生じた場合, プログラムは強制終了します. 
                              ! 引数 *err* が与えられる場合, 
                              ! プログラムは強制終了せず, 代わりに
                              ! *err* に .true. が代入されます. 
                              !
                              ! Exception handling flag. 
                              ! By default, when error occur in 
                              ! this procedure, the program aborts. 
                              ! If this *err* argument is given, 
                              ! .true. is substituted to *err* and 
                              ! the program does not abort. 

    !-----------------------------------
    !  ヒストリファイルへのデータ出力設定
    !  Configure the settings for history data output
    character(STRING):: name = ''
                              ! 変数名. Variable identifier
    character(STRING):: longname = ''
                              ! 変数の記述的名称. Descriptive name of variables
    character(STRING), allocatable:: dims(:)
                              ! 座標軸の名称. Name of axes
    character(STRING):: units = ''
                              ! 単位. Units
    type(GT_HISTORY), pointer:: gthist =>null()
                              ! gt4_history モジュール用構造体. 
                              ! Derived type for "gt4_history" module
    character(TOKEN):: precision
                              ! ヒストリデータの精度. 
                              ! Precision of history data
    logical:: average
                              ! 出力データの平均化フラグ. 
                              ! Flag for average of output data
    type(HASH):: registered_varnames
                              ! このモジュールから出力される変数名のリスト.
                              ! 
                              ! List of names of variables output 
                              ! from this module. 
    logical:: hash_end
                              ! HASH オブジェクトの終了フラグ. 
                              ! 
                              ! End flag of "HASH" object
    logical:: allvar_invalid
                              ! 無効な変数名のチェックフラグ. 
                              ! 
                              ! Check flag of invalid variable names
    character(STRING):: names_invalid
                              ! 無効な変数名. 
                              ! 
                              ! Invalid variable names
    character(STRING):: nmlfile_msg
                              ! NAMELIST ファイル名に関するメッセージ. 
                              ! 
                              ! Messages about NAMELIST file name

    !-----------------------------------
    !  作業変数
    !  Work variables
    real(DP):: UpperPressCheck
                              ! 冷却範囲の上端の気圧. 
                              ! Pressure in the top within the range of cooling
    real(DP):: LowerPressCheck
                              ! 冷却範囲の下端の気圧. 
                              ! Pressure at bottom within the range of cooling

    real(DP), parameter:: PI = 3.1415926535897930_DP
                                ! $ \pi $ . 円周率. Circular constant

    integer:: stat
    character(STRING):: cause_c
    character(*), parameter:: subname = 'PhySpoLayCreate'
  continue
    call BeginSub( subname, version )
    stat = DC_NOERR
    cause_c = ''

    !-----------------------------------------------------------------
    !  初期設定のチェック
    !  Check initialization
    !-----------------------------------------------------------------
    if ( phy_spo_lay % initialized ) then
      stat = DC_EALREADYINIT
      cause_c = 'PHYSPOLAY'
      goto 999
    end if

    !-----------------------------------------------------------------
    !  引数の正当性のチェック
    !  Validate arguments
    !-----------------------------------------------------------------
    if (imax < 1) then
      stat = DC_ENEGATIVE
      cause_c = 'imax'
      goto 999
    end if
    if (jmax < 1) then
      stat = DC_ENEGATIVE
      cause_c = 'jmax'
      goto 999
    end if
    if (kmax < 1) then
      stat = DC_ENEGATIVE
      cause_c = 'kmax'
      goto 999
    end if
    if (DelTime < 0.0_DP) then
      stat = DC_ENEGATIVE
      cause_c = 'DelTime'
      goto 999
    end if
    if ( present(DampingTime) ) then
      if (DampingTime < 0.0_DP) then
        stat = DC_ENEGATIVE
        cause_c = 'DampingTime'
        goto 999
      end if
    end if
    if ( present(UpperPress) ) then
      if (UpperPress < 0.0_DP) then
        stat = DC_ENEGATIVE
        cause_c = 'UpperPress'
        goto 999
      end if
      UpperPressCheck = UpperPress
    else
      UpperPressCheck = phy_spo_lay % UpperPress
    end if
    if ( present(LowerPress) ) then
      if (LowerPress < 0.0_DP) then
        stat = DC_ENEGATIVE
        cause_c = 'LowerPress'
        goto 999
      end if
      LowerPressCheck = LowerPress
    else
      LowerPressCheck = phy_spo_lay % LowerPress
    end if
    if ( UpperPressCheck > LowerPressCheck ) then
      call MessageNotify( 'W', subname, 'UpperPress=<%f> must be smaller than LowerPress=<%f>.', d = (/ UpperPressCheck, LowerPressCheck /) )
      stat = DCPAM_ESMALLVAL
      cause_c = 'UpperPress'
      goto 999
    end if

    !-----------------------------------------------------------------
    !  波数・格子点の設定
    !  Configure wave number and grid point
    !-----------------------------------------------------------------
    phy_spo_lay % imax  = imax 
    phy_spo_lay % jmax  = jmax 
    phy_spo_lay % kmax  = kmax 

    !-----------------------------------------------------------------
    !  座標軸の設定
    !  Configure axes
    !-----------------------------------------------------------------
    allocate( phy_spo_lay % x_Lon (0:imax-1) )
    phy_spo_lay % x_Lon = x_Lon

    allocate( phy_spo_lay % y_Lat (0:jmax-1) )
    phy_spo_lay % y_Lat = y_Lat

    allocate( phy_spo_lay % z_Sigma (0:kmax-1) )
    phy_spo_lay % z_Sigma = z_Sigma

    allocate( phy_spo_lay % x_Lon_Weight (0:imax-1) )
    if ( present(x_Lon_Weight) ) then
      phy_spo_lay % x_Lon_Weight = x_Lon_Weight
    else
      phy_spo_lay % x_Lon_Weight = 2.0_DP * PI / imax
    end if

    allocate( phy_spo_lay % y_Lat_Weight (0:jmax-1) )
    if ( present(y_Lat_Weight) ) then
      phy_spo_lay % y_Lat_Weight = y_Lat_Weight
    else
      phy_spo_lay % y_Lat_Weight = 2.0_DP / jmax
    end if

    allocate( phy_spo_lay % z_DelSigma (0:kmax-1) )
    if ( present(z_DelSigma) ) then
      phy_spo_lay % z_DelSigma = z_DelSigma
    else
      phy_spo_lay % z_DelSigma = 1.0_DP / kmax
    end if

    !-----------------------------------------------------------------
    !  その他のパラメータの設定
    !  Other parameters are configured
    !-----------------------------------------------------------------
    if ( present(DampingTime) ) then
      phy_spo_lay % DampingTime = DampingTime
    end if

    if ( present(UpperPress) ) then
      phy_spo_lay % UpperPress = UpperPress
    end if

    if ( present(LowerPress) ) then
      phy_spo_lay % LowerPress = LowerPress
    end if

    !-----------------------------------------------------------------
    !  "intavr_operate" のオブジェクトの初期設定
    !  Initialization of object of "intavr_operate"
    !-----------------------------------------------------------------
    call IntAvrOprCreate( intavr_opr = phy_spo_lay % intavr_opr, imax = imax, jmax = jmax, kmax = kmax, PI = PI, x_Lon_Weight = phy_spo_lay % x_Lon_Weight, y_Lat_Weight = phy_spo_lay % y_Lat_Weight, err = err )                                       ! (out)

    !-----------------------------------------------------------------
    !  時刻管理
    !  Time control
    !-----------------------------------------------------------------
    if ( present(current_time_value) .and. present(current_time_unit) ) then
      call DCDiffTimeCreate( diff = phy_spo_lay % current_time, value = real( current_time_value, DP ), unit = current_time_unit )                ! (in)
    else
      call DCDiffTimeCreate( diff = phy_spo_lay % current_time, value = 0.0_DP, unit = 'sec' )       ! (in)
    end if

    call DCDiffTimeCreate( diff = phy_spo_lay % delta_time, value = DelTime, unit = 'sec' )                 ! (in)

    !-----------------------------------------------------------------
    !  ヒストリファイルへのデータ出力設定
    !  Configure the settings for history data output
    !-----------------------------------------------------------------

    !-------------------------
    !  デフォルト値
    !  Default values
    call HstNmlInfoCreate( gthstnml = phy_spo_lay % gthstnml ) ! (inout)

    !-------------------------
    !  オプショナル引数からの値
    !  Values from optional arguments
    call HstNmlInfoAdd( gthstnml = phy_spo_lay % gthstnml, name = '', interval_value = history_interval_value, interval_unit = history_interval_unit, precision = history_precision, average = .false., fileprefix = history_fileprefix )          ! (in)

    if ( present(history_varlist) ) then
      call HstNmlInfoAdd( gthstnml = phy_spo_lay % gthstnml, name = history_varlist )                 ! (in)
    end if

    !-----------------------------------------------------------------
    !  NAMELIST からの値の読み込み
    !  Load values from NAMELIST
    !-----------------------------------------------------------------
    if ( present_and_not_empty(nmlfile) ) then
      call MessageNotify( 'M', subname, 'Loading NAMELIST file "%c" ...', c1 = trim(nmlfile) )
      call NmlRead ( nmlfile = nmlfile, DampingTime = phy_spo_lay % DampingTime, UpperPress  = phy_spo_lay % UpperPress, LowerPress  = phy_spo_lay % LowerPress, gthstnml    = phy_spo_lay % gthstnml, err = err )                              ! (out)
      if ( present_and_true(err) ) then
        call MessageNotify( 'W', subname, '"%c" can not be read.', c1 = trim(nmlfile) )
        stat = DC_ENOFILEREAD
        cause_c = nmlfile
        goto 999
      end if
    end if

    call HstNmlInfoEndDefine( gthstnml = phy_spo_lay % gthstnml, err = err )                              ! (out)
    if ( present_and_true( err ) ) then
      stat = USR_ERRNO
      goto 999
    end if


    !-----------------------------------------------------------------
    !  データ出力の初期設定
    !  Initialize data output
    !-----------------------------------------------------------------

    !-------------------------
    !  xyz_DUDtSpoDamping の出力設定
    !  Configure the settings for "xyz_DUDtSpoDamping" output
    name = 'DUDtSpoDamping'
    longname = 'zonal wind tendency by damping in sponge layer'
    units = 'm s-2'
    allocate( dims(4) )
    dims = StoA( 'lon', 'lat', 'sig', 'time' )

    ! 出力ファイルの初期設定.
    !   * gthist (gt4_history#GT_HISTORY) が設定される.
    ! Initialize output file.
    !   * "gthist" (gt4_history#GT_HISTORY) is configured.
    call output_init  ! これは内部サブルーチン. This is an internal subroutine
    deallocate( dims )

    !-------------------------
    !  xyz_DVDtSpoDamping の出力設定
    !  Configure the settings for "xyz_DVDtSpoDamping" output
    name = 'DVDtSpoDamping'
    longname = 'meridional wind tendency by damping in sponge layer'
    units = 'm s-2'
    allocate( dims(4) )
    dims = StoA( 'lon', 'lat', 'sig', 'time' )

    ! 出力ファイルの初期設定.
    !   * gthist (gt4_history#GT_HISTORY) が設定される.
    ! Initialize output file.
    !   * "gthist" (gt4_history#GT_HISTORY) is configured.
    call output_init  ! これは内部サブルーチン. This is an internal subroutine
    deallocate( dims )

!!$    !-------------------------
!!$    !  y_Data2 の出力設定
!!$    !  Configure the settings for "y_Data2" output
!!$    name = 'Data2'
!!$    longname = 'Sample data (2)'
!!$    units = '1'
!!$    allocate( dims(2) )
!!$    dims = StoA( 'lat', 'time' )
!!$
!!$    ! 出力ファイルの初期設定.
!!$    !   * gthist (gt4_history#GT_HISTORY) が設定される.
!!$    ! Initialize output file.
!!$    !   * "gthist" (gt4_history#GT_HISTORY) is configured.
!!$    call output_init  ! これは内部サブルーチン. This is an internal subroutine
!!$    ! 属性の付加などを行う場合には以下のようにする. 
!!$    ! Describe codes as follows in order to add attributes etc.
!!$    if ( associated( gthist ) ) then
!!$      call HistoryAddAttr( &
!!$        & history = gthist, &                           ! (inout)
!!$        & varname = name, attrname = 'missing_value', & ! (in)
!!$        &   value = 9999.0 )                            ! (in)
!!$    end if
!!$    deallocate( dims )

    !-----------------------------------------------------------------
    !  このモジュールから出力される変数名のリストを表示
    !  Print list of names of variables output from this module
    !-----------------------------------------------------------------
    call Printf( STDOUT, ' *** MESSAGE *** +---- "%c" module output varnames list -----', c1 = 'phy_sponge_layer' )
    call DCHashRewind( hashv = registered_varnames ) ! (inout)
    do
      call DCHashNext( hashv = registered_varnames, key = name, value = longname, end = hash_end ) ! (out)
      if ( hash_end ) exit
      call Printf( STDOUT, ' *** MESSAGE *** |      "%c" (%c)', c1 = trim(name), c2 = trim(longname) )
    enddo
    call DCHashDelete( hashv = registered_varnames ) ! (inout)
    call Printf( STDOUT, ' *** MESSAGE *** `----------------------------------------' )

    !-----------------------------------------------------------------
    !  無効な変数名のチェック
    !  Check invalid variable names
    !-----------------------------------------------------------------
    call HstNmlInfoAllNameValid( gthstnml = phy_spo_lay % gthstnml, invalid = allvar_invalid, names = names_invalid, err = err )                                        ! (out)

    if ( allvar_invalid ) then
      stat = HST_EBADVARNAME
      cause_c = names_invalid
      phy_spo_lay % initialized = .true.
      call PhySpoLayClose( phy_spo_lay, err = err )                     ! (out)
      if ( present(nmlfile) ) then
        nmlfile_msg = ' or NAMELIST "phy_sponge_layer_history_nml" in "'// trim(nmlfile) // '"'
      else
        nmlfile_msg = ''
      end if
      call MessageNotify( 'W', subname, 'names "%c" from an optional argument "history_varlist"%c are invalid.', c1 = trim(names_invalid), c2 = trim(nmlfile_msg) )
      goto 999
    end if

    !-----------------------------------------------------------------
    !  設定値の正当性のチェック
    !  Validate setting values
    !-----------------------------------------------------------------
!!$    if ( phy_spo_lay % CoefAlpha < 0.0_DP ) then
!!$      stat = DC_ENEGATIVE
!!$      cause_c = 'CoefAlpha'
!!$      goto 999
!!$    end if

    !-----------------------------------------------------------------
    !  終了処理, 例外処理
    !  Termination and Exception handling
    !-----------------------------------------------------------------
    phy_spo_lay % initialized = .true.
999 continue
    call StoreError( stat, subname, err, cause_c )
    call EndSub( subname )

  contains

    subroutine output_init
      !
      ! 変数 *name* に関して出力ファイルの初期設定を行います. 
      ! 出力ファイル名や出力間隔などの情報は phy_spo_lay % gthstnml 
      ! から取り出されます. 
      ! 
      ! 変数 *name* に関して出力が行われる場合には, 
      ! *gthist* に出力先ファイルの gt4_history#GT_HISTORY
      ! 型変数を結合させます. そうでない場合は, *gthist* を空状態にします. 
      !
      ! また, 出力データの精度を precision に, 
      ! 出力データ平均化の可否を average に設定します. 
      ! 
      ! 標準出力に表示される変数リスト *registered_varnames* に
      ! *name*, *longname*, *dims*, *units* が登録されます. 
      !
      ! An output file is initialized for a variable *name*. 
      ! Information such as the output filename and output intervals 
      ! is taken out of "phy_spo_lay % gthstnml". 
      !
      ! When output is done for the variable *name*, *gthist* is 
      ! associated with the "gt4_history#GT_HISTORY" variable of 
      ! the output file. Otherwise, *gthist* is nullified. 
      !
      ! Moreover, the accuracy of output data is set to *precision*, and 
      ! right or wrong of averaging the output data is set to *average*. 
      !
      ! *name*, *longname*, *dims*, *units* are registered to 
      ! a list of variables *registered_varnames* that is printed to 
      ! standard output. 
      !
      use dc_date, only: DCDiffTimeCreate, EvalSec, EvalByUnit
      use gt4_history_nmlinfo, only: GTHST_NMLINFO, HstNmlInfoInitialized, HstNmlInfoInquire, HstNmlInfoOutputValid, HstNmlInfoAssocGtHist, HstNmlInfoPutLine, HstNmlInfoSetValidName
      use gt4_history, only: GT_HISTORY, HistoryCreate, HistoryAddVariable, HistoryPut, HistoryAddAttr, HistoryInitialized

      !-----------------------------------
      !  作業変数
      !  Work variables
      character(STRING):: file
                                ! ヒストリデータのファイル名. 
                                ! History data filenames
      character(STRING):: dims_str
                                ! 座標軸のリスト. 
                                ! List of axes
      real:: interval_value
                                ! ヒストリデータの出力間隔の数値. 
                                ! Numerical value for interval of history data output
      character(TOKEN):: interval_unit
                                ! ヒストリデータの出力間隔の単位. 
                                ! Unit for interval of history data output
    continue
      !-----------------------------------------------------------------
      !  標準出力に表示される変数の登録
      !  Register a variable name for print to standard output
      !-----------------------------------------------------------------
      if ( allocated(dims) ) then
        dims_str = JoinChar( dims, ',' )
      else
        dims_str = ''
      end if
      call DCHashPut( hashv = registered_varnames, key = name, value = trim( longname ) // ' [' // trim( units ) // '] {' // trim( dims_str ) // '}' ) ! (in)

      !-----------------------------------------------------------------
      !  変数の初期化
      !  Initialize variable
      !-----------------------------------------------------------------
      nullify( gthist )
      precision = 'float'
      average = .false.

      !-----------------------------------------------------------------
      !  変数名の有効性を設定
      !  Set validation of the variable name
      !-----------------------------------------------------------------
      call HstNmlInfoSetValidName( gthstnml = phy_spo_lay % gthstnml, name = name )                        ! (in)

      !-----------------------------------------------------------------
      !  出力が有効かどうかを確認する
      !  Confirm whether the output is effective
      !-----------------------------------------------------------------
      if ( .not. HstNmlInfoOutputValid( phy_spo_lay % gthstnml, name ) ) then
        return
      end if

      !-----------------------------------------------------------------
      !  GT_HISTORY 変数の取得
      !  Get "GT_HISTORY" variable
      !-----------------------------------------------------------------
      call HstNmlInfoAssocGtHist( gthstnml = phy_spo_lay % gthstnml, name = name, history = gthist, err = err )                          ! (out)
      if ( present_and_true( err ) ) return

      call HstNmlInfoInquire( gthstnml = phy_spo_lay % gthstnml, name = name, precision = precision, average = average, err = err )                          ! (out)
      if ( present_and_true( err ) ) return

      !-----------------------------------------------------------------
      !  GT_HISTORY 変数の初期設定の確認
      !  Check initialization of "GT_HISTORY" variable
      !-----------------------------------------------------------------
      if ( HistoryInitialized( gthist ) ) then

        !---------------------------------------------------------------
        !  HistoryAddVariable による変数作成
        !  A variable is created by "HistoryAddVariable"
        !---------------------------------------------------------------
        call HistoryAddVariable( history = gthist, varname = name,         dims = dims, longname = longname,    units = units, xtype = precision, average = average )  ! (in)
        return
      end if

      !-----------------------------------------------------------------
      !  HistoryCreate のための設定値の取得
      !  Get the settings for "HistoryCreate"
      !-----------------------------------------------------------------
      call HstNmlInfoInquire( gthstnml = phy_spo_lay % gthstnml, name = name, file = file, interval_unit = interval_unit, interval_value = interval_value, err = err )                          ! (out)
      if ( present_and_true( err ) ) return

      !-----------------------------------------------------------------
      !  HistoryCreate によるファイル作成
      !  Files are created by "HistoryCreate"
      !-----------------------------------------------------------------
      call HistoryCreate( history = gthist, file = file, title = 'Sponge layer for damping velocity', source = 'dcpam4 project : ' // trim(version), institution = 'GFD Dennou Club', dims = StoA( 'lon', 'lat', 'sig', 'time' ), dimsizes = (/ phy_spo_lay % imax, phy_spo_lay % jmax, phy_spo_lay % kmax, 0 /), longnames = StoA('longitude', 'latitude', 'sigma at layer midpoints', 'time'), units = StoA( 'degree_east', 'degree_north', '1', interval_unit ), origin = real( EvalbyUnit( phy_spo_lay % current_time, interval_unit) ), interval = interval_value, err = err )                                            ! (out)
      if ( present_and_true( err ) ) then
        nullify( gthist )
        return
      end if

      call HistoryAddAttr( history = gthist, varname = 'lon', attrname = 'standard_name', value = 'longitude' )                          ! (in)
      call HistoryAddAttr( history = gthist, varname = 'lat', attrname = 'standard_name', value = 'latitude' )                           ! (in)
      call HistoryAddAttr( history = gthist, varname = 'sig', attrname = 'standard_name', value = 'atmosphere_sigma_coordinate' )          ! (in)
      call HistoryAddAttr( history = gthist, varname = 'time', attrname = 'standard_name', value = 'time' )                                ! (in)
      call HistoryAddAttr( history = gthist, varname = 'sig', attrname = 'positive', value = 'down' )                                 ! (in)

      call HistoryPut( history = gthist, varname = 'lon', array = phy_spo_lay % x_Lon / PI * 180.0_DP ) ! (in)
      call HistoryPut( history = gthist, varname = 'lat', array = phy_spo_lay % y_Lat / PI * 180.0_DP ) ! (in)
      call HistoryPut( history = gthist, varname = 'sig', array = phy_spo_lay % z_Sigma ) ! (in)

      call HistoryAddVariable( history = gthist, varname = 'lon_weight', dims = StoA('lon'), longname = 'weight for integration in longitude', units = 'radian', xtype = 'double' )                ! (in)
      call HistoryAddAttr( history = gthist, varname = 'lon', attrname = 'gt_calc_weight', value = 'lon_weight' )                              ! (in)
      call HistoryPut( history = gthist, varname = 'lon_weight', array = phy_spo_lay % x_Lon_Weight )    ! (in)

      call HistoryAddVariable( history = gthist, varname = 'lat_weight', dims = StoA('lat'), longname = 'weight for integration in latitude', units = 'radian', xtype = 'double' )                ! (in)
      call HistoryAddAttr( history = gthist, varname = 'lat', attrname = 'gt_calc_weight', value = 'lat_weight' )                              ! (in)
      call HistoryPut( history = gthist, varname = 'lat_weight', array = phy_spo_lay % y_Lat_Weight )    ! (in)

      call HistoryAddVariable( history = gthist, varname = 'sig_weight', dims = StoA('sig'), longname = 'weight for integration in sigma', units = '1', xtype = 'float' )                      ! (in)
      call HistoryAddAttr( history = gthist, varname = 'sig', attrname = 'gt_calc_weight', value = 'sig_weight' )                              ! (in)
      call HistoryPut( history = gthist, varname = 'sig_weight', array = phy_spo_lay % z_DelSigma )      ! (in)

      !-----------------------------------------------------------------
      !  HistoryAddVariable による変数作成
      !  A variable is created by "HistoryAddVariable"
      !-----------------------------------------------------------------
      if ( HistoryInitialized( gthist ) ) then
        call HistoryAddVariable( history = gthist, varname = name,         dims = dims, longname = longname,    units = units, xtype = precision, average = average )  ! (in)
      else
        nullify( gthist )
      end if

    end subroutine output_init

  end subroutine PhySpoLayCreate
Function :
result :logical
phy_spo_lay :type(PHYSPOLAY), intent(in)

phy_spo_lay が初期設定されている場合には .true. が, 初期設定されていない場合には .false. が返ります.

If phy_spo_lay is initialized, .true. is returned. If phy_spo_lay is not initialized, .false. is returned.

[Source]

  logical function PhySpoLayInitialized( phy_spo_lay ) result(result)
    !
    ! *phy_spo_lay* が初期設定されている場合には .true. が, 
    ! 初期設定されていない場合には .false. が返ります. 
    !
    ! If *phy_spo_lay* is initialized, .true. is returned. 
    ! If *phy_spo_lay* is not initialized, .false. is returned. 
    !
    implicit none
    type(PHYSPOLAY), intent(in):: phy_spo_lay
  continue
    result = phy_spo_lay % initialized
  end function PhySpoLayInitialized
Subroutine :
phy_spo_lay :type(PHYSPOLAY), intent(in)
unit :integer, intent(in), optional
: 出力先の装置番号. デフォルトの出力先は標準出力.

Unit number for output. Default value is standard output.

indent :character(*), intent(in), optional
: 表示されるメッセージの字下げ.

Indent of displayed messages.

err :logical, intent(out), optional
: 例外処理用フラグ. デフォルトでは, この手続き内でエラーが 生じた場合, プログラムは強制終了します. 引数 err が与えられる場合, プログラムは強制終了せず, 代わりに err に .true. が代入されます.

Exception handling flag. By default, when error occur in this procedure, the program aborts. If this err argument is given, .true. is substituted to err and the program does not abort.

引数 phy_spo_lay に設定されている情報を印字します. デフォルトではメッセージは標準出力に出力されます. unit に装置番号を指定することで, 出力先を変更することが可能です.

Print information of phy_spo_lay. By default messages are output to standard output. Unit number for output can be changed by unit argument.

[Source]

  subroutine PhySpoLayPutLine( phy_spo_lay, unit, indent, err )
    !
    ! 引数 *phy_spo_lay* に設定されている情報を印字します. 
    ! デフォルトではメッセージは標準出力に出力されます. 
    ! *unit* に装置番号を指定することで, 出力先を変更することが可能です. 
    !
    ! Print information of *phy_spo_lay*. 
    ! By default messages are output to standard output. 
    ! Unit number for output can be changed by *unit* argument. 
    !
    use dc_trace, only: BeginSub, EndSub
    use dc_string, only: PutLine, Printf, Split, StrInclude, StoA, JoinChar
    use dc_types, only: DP, STRING, TOKEN, STDOUT
    use dc_date, only: EvalSec
    use dc_error, only: StoreError, DC_NOERR, DC_ENOTINIT
    use gt4_history_nmlinfo, only: GTHST_NMLINFO, HstNmlInfoPutLine
    use intavr_operate, only: IntAvrOprPutLine
    implicit none
    type(PHYSPOLAY), intent(in):: phy_spo_lay
    integer, intent(in), optional:: unit
                              ! 出力先の装置番号. 
                              ! デフォルトの出力先は標準出力. 
                              !
                              ! Unit number for output. 
                              ! Default value is standard output. 
    character(*), intent(in), optional:: indent
                              ! 表示されるメッセージの字下げ. 
                              !
                              ! Indent of displayed messages. 
    logical, intent(out), optional:: err
                              ! 例外処理用フラグ. 
                              ! デフォルトでは, この手続き内でエラーが
                              ! 生じた場合, プログラムは強制終了します. 
                              ! 引数 *err* が与えられる場合, 
                              ! プログラムは強制終了せず, 代わりに
                              ! *err* に .true. が代入されます. 
                              !
                              ! Exception handling flag. 
                              ! By default, when error occur in 
                              ! this procedure, the program aborts. 
                              ! If this *err* argument is given, 
                              ! .true. is substituted to *err* and 
                              ! the program does not abort. 

    !-----------------------------------
    !  作業変数
    !  Work variables
    integer:: stat
    character(STRING):: cause_c
    integer:: out_unit
    integer:: indent_len
    character(STRING):: indent_str
    character(*), parameter:: subname = 'PhySpoLayPutLine'
  continue
    call BeginSub( subname )
    stat = DC_NOERR
    cause_c = ''

    !-----------------------------------------------------------------
    !  出力先装置番号と字下げの設定
    !  Configure output unit number and indents
    !-----------------------------------------------------------------
    if ( present(unit) ) then
      out_unit = unit
    else
      out_unit = STDOUT
    end if

    indent_len = 0
    indent_str = ''
    if ( present(indent) ) then
      if ( len(indent) /= 0 ) then
        indent_len = len(indent)
        indent_str(1:indent_len) = indent
      end if
    end if

    !-----------------------------------------------------------------
    !  "PHYSPOLAY" の設定の印字
    !  Print the settings for "PHYSPOLAY"
    !-----------------------------------------------------------------
    if ( phy_spo_lay % initialized ) then
      call Printf( unit = out_unit, fmt = indent_str(1:indent_len) // '#<PHYSPOLAY:: @initialized=%y', l = (/ phy_spo_lay % initialized /) )   ! (in)

      call Printf( unit = out_unit, fmt = indent_str(1:indent_len) // ' @imax=%d @jmax=%d @kmax=%d', i = (/ phy_spo_lay % imax, phy_spo_lay % jmax, phy_spo_lay % kmax /) )          ! (in)

      call Printf( unit = out_unit, fmt = indent_str(1:indent_len) // ' @current_time=%f [sec] @delta_time=%f [sec]', d = (/ EvalSec( phy_spo_lay % current_time ), EvalSec( phy_spo_lay % delta_time ) /) )      ! (in)

      call Printf( unit = out_unit, fmt = indent_str(1:indent_len) // ' @DampingTime=%f', d = (/phy_spo_lay % DampingTime/) ) ! (in)

      call Printf( unit = out_unit, fmt = indent_str(1:indent_len) // ' @UpperPress=%f, @LowerPress=%f', d = (/phy_spo_lay % UpperPress, phy_spo_lay % LowerPress/) ) ! (in)

      call PutLine( array = phy_spo_lay % x_Lon, unit = out_unit, lbounds = lbound(phy_spo_lay % x_Lon), ubounds = ubound(phy_spo_lay % x_Lon), indent = indent_str(1:indent_len) // ' @x_Lon=' )                        ! (in)

      call PutLine( array = phy_spo_lay % y_Lat, unit = out_unit, lbounds = lbound(phy_spo_lay % y_Lat), ubounds = ubound(phy_spo_lay % y_Lat), indent = indent_str(1:indent_len) // ' @y_Lat=' )                        ! (in)

      call PutLine( array = phy_spo_lay % z_Sigma, unit = out_unit, lbounds = lbound(phy_spo_lay % z_Sigma), ubounds = ubound(phy_spo_lay % z_Sigma), indent = indent_str(1:indent_len) // ' @z_Sigma=' )                    ! (in)

      call IntAvrOprPutLine( phy_spo_lay % intavr_opr, unit = out_unit, indent = indent_str(1:indent_len) // ' ', err = err )                                    ! (out)

      call Printf( unit = out_unit, fmt = indent_str(1:indent_len) // ' @gthstnml=' )               ! (in)
      call HstNmlInfoPutLine( gthstnml = phy_spo_lay % gthstnml, unit = out_unit, indent = indent_str(1:indent_len) // '  ' ) ! (in)

      call Printf( unit = out_unit, fmt = indent_str(1:indent_len) // '>' ) ! (in)
    else
      call Printf( unit = out_unit, fmt = indent_str(1:indent_len) // '#<PHYSPOLAY:: @initialized=%y>', l = (/phy_spo_lay % initialized/) )       ! (in)
    end if

    !-----------------------------------------------------------------
    !  終了処理, 例外処理
    !  Termination and Exception handling
    !-----------------------------------------------------------------
999 continue
    call StoreError( stat, subname, err, cause_c )
    call EndSub( subname )

  end subroutine PhySpoLayPutLine
Subroutine :
phy_spo_lay :type(PHYSPOLAY), intent(inout)
current_time_value :real, intent(in)
: 現在時刻の数値. Numerical value of current time
current_time_unit :character(*), intent(in)
: 現在時刻の単位. Unit of current time
err :logical, intent(out), optional
: 例外処理用フラグ. デフォルトでは, この手続き内でエラーが 生じた場合, プログラムは強制終了します. 引数 err が与えられる場合, プログラムは強制終了せず, 代わりに err に .true. が代入されます.

Exception handling flag. By default, when error occur in this procedure, the program aborts. If this err argument is given, .true. is substituted to err and the program does not abort.

phy_spo_lay に対して時刻の設定を行います.

ヒストリデータを出力している場合には, ヒストリデータの出力時刻も設定します. 一度でもこのサブルーチンを呼んだ場合には, それ以後のヒストリデータ出力前にこのサブルーチンを呼び出し, 時刻の設定を行ってください. また, データを出力するサブルーチンに対しても オプショナル引数 historyput_flag に .true. を与えてください.

なお, 与えられた phy_spo_layPhySpoLayCreate によって初期設定されていない場合, プログラムはエラーを発生させます.

Set time to phy_spo_lay.

When history data are output, the output time of history data are specified. Once this subroutine is called, the time of history data must be specified by this routine before history data output. In additional, give ".true." to an optional argument "historyput_flag" of a data output subroutine.

If phy_spo_lay is not initialized by "PhySpoLayCreate" yet, error is occurred.

[Source]

  subroutine PhySpoLaySetTime( phy_spo_lay, current_time_value, current_time_unit, err )
    !
    ! *phy_spo_lay* に対して時刻の設定を行います. 
    !
    ! ヒストリデータを出力している場合には, 
    ! ヒストリデータの出力時刻も設定します. 
    ! 一度でもこのサブルーチンを呼んだ場合には, 
    ! それ以後のヒストリデータ出力前にこのサブルーチンを呼び出し, 
    ! 時刻の設定を行ってください. 
    ! また, データを出力するサブルーチンに対しても 
    ! オプショナル引数 historyput_flag に .true. を与えてください. 
    !
    ! なお, 与えられた *phy_spo_lay* が PhySpoLayCreate 
    ! によって初期設定されていない場合, プログラムはエラーを発生させます. 
    !
    ! Set time to *phy_spo_lay*. 
    !
    ! When history data are output, 
    ! the output time of history data are specified. 
    ! Once this subroutine is called, the time of history data must be 
    ! specified by this routine before history data output. 
    ! In additional, give ".true." to an optional argument 
    ! "historyput_flag" of a data output subroutine. 
    !
    ! If *phy_spo_lay* is not initialized 
    ! by "PhySpoLayCreate" yet, error is occurred. 
    !
    use dc_trace, only: BeginSub, EndSub
    use dc_string, only: PutLine, Printf, Split, StrInclude, StoA, JoinChar
    use dc_types, only: DP, STRING, TOKEN, STDOUT
    use dc_date, only: DCDiffTimeCreate, EvalbyUnit
    use dc_error, only: StoreError, DC_NOERR, DC_ENOTINIT
    use gt4_history_nmlinfo, only: HstNmlInfoAdd, HstNmlInfoInquire, HstNmlInfoNames, HstNmlInfoAssocGtHist, HstNmlInfoOutputStepDisable, HstNmlInfoPutLine
    use gt4_history, only: GT_HISTORY, HistorySetTime, HistoryInitialized
    implicit none
    type(PHYSPOLAY), intent(inout):: phy_spo_lay
    real, intent(in):: current_time_value
                              ! 現在時刻の数値. Numerical value of current time
    character(*), intent(in):: current_time_unit
                              ! 現在時刻の単位. Unit of current time
    logical, intent(out), optional:: err
                              ! 例外処理用フラグ. 
                              ! デフォルトでは, この手続き内でエラーが
                              ! 生じた場合, プログラムは強制終了します. 
                              ! 引数 *err* が与えられる場合, 
                              ! プログラムは強制終了せず, 代わりに
                              ! *err* に .true. が代入されます. 
                              !
                              ! Exception handling flag. 
                              ! By default, when error occur in 
                              ! this procedure, the program aborts. 
                              ! If this *err* argument is given, 
                              ! .true. is substituted to *err* and 
                              ! the program does not abort. 

    !-----------------------------------
    !  ヒストリファイルへのデータ出力設定
    !  Configure the settings for history data output
    character(STRING):: name = ''
                              ! 変数名. Variable identifier
    character(TOKEN):: interval_unit
                              ! ヒストリデータの出力間隔の単位. 
                              ! Unit for interval of history data output
    character(STRING):: varnames
                              ! 変数名リスト. 
                              ! List of variables
    character(TOKEN), pointer:: varnames_array(:) =>null()
                              ! 変数名リスト配列. 
                              ! List of variables (array) 
    integer:: i, vnmax
    type(GT_HISTORY), pointer:: gthist =>null()
                              ! gt4_history モジュール用構造体. 
                              ! Derived type for "gt4_history" module

    !-----------------------------------
    !  作業変数
    !  Work variables
    integer:: stat
    character(STRING):: cause_c
    character(*), parameter:: subname = 'PhySpoLaySetTime'
  continue
    call BeginSub( subname )
    stat = DC_NOERR
    cause_c = ''

    !-----------------------------------------------------------------
    !  初期設定のチェック
    !  Check initialization
    !-----------------------------------------------------------------
    if ( .not. phy_spo_lay % initialized ) then
      stat = DC_ENOTINIT
      cause_c = 'PHYSPOLAY'
      goto 999
    end if

    !-----------------------------------------------------------------
    !  時刻設定
    !  Configure time
    !-----------------------------------------------------------------
    call DCDiffTimeCreate( diff = phy_spo_lay % current_time, value = real( current_time_value, DP ), unit = current_time_unit )                ! (in)

    !-----------------------------------------------------------------
    !  ヒストリファイルへのデータの時刻設定
    !  Configure the time of history data
    !-----------------------------------------------------------------
    varnames = HstNmlInfoNames( phy_spo_lay % gthstnml )
    call Split( str = varnames, sep = ',', carray = varnames_array )            ! (out)
    vnmax = size( varnames_array )

    do i = 1, vnmax
      name = varnames_array(i)
      if ( trim( name ) == '' ) exit

      call HstNmlInfoOutputStepDisable( gthstnml = phy_spo_lay % gthstnml, name = name, err = err )                              ! (out)

      nullify( gthist )
      call HstNmlInfoAssocGtHist( gthstnml = phy_spo_lay % gthstnml, name = name, history = gthist, err = err )                              ! (out)

      if ( HistoryInitialized( gthist ) ) then
        call HstNmlInfoInquire( gthstnml = phy_spo_lay % gthstnml, name = name, interval_unit = interval_unit )          ! (out)

        call HistorySetTime( history = gthist, time = real( EvalbyUnit( phy_spo_lay % current_time, interval_unit) ) ) ! (in)
      end if

    end do

    !-----------------------------------------------------------------
    !  終了処理, 例外処理
    !  Termination and Exception handling
    !-----------------------------------------------------------------
999 continue
    call StoreError( stat, subname, err, cause_c )
    call EndSub( subname )

  end subroutine PhySpoLaySetTime
PhySpoLayinitialized( phy_spo_lay ) result(result)
Function :
result :logical
phy_spo_lay :type(PHYSPOLAY), intent(in)

phy_spo_lay が初期設定されている場合には .true. が, 初期設定されていない場合には .false. が返ります.

If phy_spo_lay is initialized, .true. is returned. If phy_spo_lay is not initialized, .false. is returned.

Alias for PhySpoLayInitialized

Private Instance methods

NmlRead( nmlfile, DampingTime, UpperPress, LowerPress, gthstnml, [err] )
Subroutine :
nmlfile :character(*), intent(in)
: NAMELIST ファイルの名称. NAMELIST file name
DampingTime :real(DP), intent(inout)
: 減衰の時定数 (単位: 秒). Damping time constant (units: seconds)
UpperPress :real(DP), intent(inout)
: 減衰範囲の上端の気圧. Pressure in the top within the range of damping
LowerPress :real(DP), intent(inout)
: 減衰範囲の下端の気圧. Pressure at bottom within the range of damping
gthstnml :type(GTHST_NMLINFO), intent(inout)
: NAMELIST#phy_sponge_layer_history_nml から入手される個別のデータ出力情報.

初期設定やデフォルト値の設定などを 行った後に与えること.

Individual data output information from "NAMELIST#phy_sponge_layer_history_nml".

Before this argument is given to this procedure, initialize and configure the defaut settings.

err :logical, intent(out), optional
: 例外処理用フラグ. デフォルトでは, この手続き内でエラーが 生じた場合, プログラムは強制終了します. 引数 err が与えられる場合, プログラムは強制終了せず, 代わりに err に .true. が代入されます.

Exception handling flag. By default, when error occur in this procedure, the program aborts. If this err argument is given, .true. is substituted to err and the program does not abort.

NAMELIST ファイル nmlfile から値を入力するための サブルーチンです. PhySpoLayCreate 内で呼び出されることを想定しています.

値が NAMELIST ファイル内で指定されていない場合には, 入力された値がそのまま返ります.

なお, nmlfile に空文字が与えられた場合, または 与えられた nmlfile を読み込むことができない場合, プログラムはエラーを発生させます.

This is a subroutine to input values from NAMELIST file nmlfile. This subroutine is expected to be called by "PhySpoLayCreate".

A value not specified in NAMELIST file is returned without change.

If nmlfile is empty, or nmlfile can not be read, error is occurred.

This procedure input/output NAMELIST#phy_sponge_layer_nml, NAMELIST#phy_sponge_layer_history_nml .

Alias for PhySpoLayNmlRead

Subroutine :
phy_spo_lay :type(PHYSPOLAY), intent(inout)
xyz_U(0:phy_spo_lay%imax-1, 0:phy_spo_lay%jmax-1, 0:phy_spo_lay%kmax-1) :real(DP), intent(in)
: $ u $ . 東西風速. Zonal wind
xyz_V(0:phy_spo_lay%imax-1, 0:phy_spo_lay%jmax-1, 0:phy_spo_lay%kmax-1) :real(DP), intent(in)
: $ v $ . 南北風速. Meridional wind
xy_Ps(0:phy_spo_lay%imax-1, 0:phy_spo_lay%jmax-1) :real(DP), intent(in)
: $ p_s $ . 地表面気圧. Surface pressure
xyz_DUDtSpoDamping(0:phy_spo_lay%imax-1, 0:phy_spo_lay%jmax-1, 0:phy_spo_lay%kmax-1) :real(DP), intent(out)
: $ DP{u}{t} $ . 減衰の効果による東西風速変化. Zonal wind tendency by damping effect
xyz_DVDtSpoDamping(0:phy_spo_lay%imax-1, 0:phy_spo_lay%jmax-1, 0:phy_spo_lay%kmax-1) :real(DP), intent(out)
: $ DP{v}{t} $ . 減衰の効果による南北風速変化. Meridional wind tendency by damping effect
historyput_flag :logical, intent(in), optional
: データ出力のフラグ. SetTime によって時刻を明示的に 指定した場合には, この引数に .true. または .false. を指定する ことでデータ出力のオンオフを 明示的に指定する必要があります. デフォルトは .false. です.

Data output flag. When time is specified by "SetTime", explicit specification of data output on/off by specifying ".true." or ".false." to this argument. Default value is ".false.".

err :logical, intent(out), optional
: 例外処理用フラグ. デフォルトでは, この手続き内でエラーが 生じた場合, プログラムは強制終了します. 引数 err が与えられる場合, プログラムは強制終了せず, 代わりに err に .true. が代入されます.

Exception handling flag. By default, when error occur in this procedure, the program aborts. If this err argument is given, .true. is substituted to err and the program does not abort.

東西風速 xyz_U, 南北風速 xyz_V および 地表面気圧 xy_Ps から, 減衰の効果による東西風速と南北風速の変化率 xyz_DUDtSpoDamping, xyz_DVDtSpoDamping を算出して返します. 計算後, phy_spo_lay 内の時刻を更新します.

なお, 与えられた phy_spo_layPhySpoLayCreate によって初期設定されていない場合, プログラムはエラーを発生させます.

Calculate and return zonal wind and meridional tendency by damping in sponge layer xyz_DUDtSpoDamping, xyz_DVDtSpoDamping from zonal wind xyz_U and meridional wind xyz_V and surface pressure xy_Ps. After the calculation, time stored in phy_spo_lay is updated.

If phy_spo_lay is not initialized by "PhySpoLayCreate" yet, error is occurred.

[Source]

  subroutine PhySpoLayDamping( phy_spo_lay, xyz_U, xyz_V, xy_Ps, xyz_DUDtSpoDamping, xyz_DVDtSpoDamping, historyput_flag, err )
    !
    ! 東西風速 *xyz_U*, 南北風速 *xyz_V* および 
    ! 地表面気圧 *xy_Ps* から, 
    ! 減衰の効果による東西風速と南北風速の変化率
    ! *xyz_DUDtSpoDamping*, *xyz_DVDtSpoDamping* 
    ! を算出して返します. 
    ! 計算後, *phy_spo_lay* 内の時刻を更新します. 
    !
    ! なお, 与えられた *phy_spo_lay* が PhySpoLayCreate 
    ! によって初期設定されていない場合, プログラムはエラーを発生させます. 
    !
    ! Calculate and return
    ! zonal wind and meridional tendency by damping in sponge layer
    ! *xyz_DUDtSpoDamping*, *xyz_DVDtSpoDamping* 
    ! from
    ! zonal wind *xyz_U* and meridional wind *xyz_V* and
    ! surface pressure *xy_Ps*.
    ! After the calculation, time stored in *phy_spo_lay* is updated. 
    !
    ! If *phy_spo_lay* is not initialized 
    ! by "PhySpoLayCreate" yet, error is occurred. 
    !
    use dc_trace, only: BeginSub, EndSub
    use dc_string, only: PutLine, Printf, Split, StrInclude, StoA, JoinChar
    use dc_types, only: DP, STRING, TOKEN, STDOUT
    use dc_present, only: present_and_true
    use dc_date, only: mod, operator(+), operator(==), EvalbyUnit, EvalSec
    use dc_error, only: StoreError, DC_NOERR, DC_ENOTINIT, USR_ERRNO
    use gt4_history_nmlinfo, only: HstNmlInfoInquire, HstNmlInfoOutputValid, HstNmlInfoAssocGtHist, HstNmlInfoPutLine
    use gt4_history, only: GT_HISTORY, HistoryPut, HistoryInitialized
    use intavr_operate, only: ya_AvrLon_xya
    implicit none
    type(PHYSPOLAY), intent(inout):: phy_spo_lay
    real(DP), intent(in):: xyz_U (0:phy_spo_lay%imax-1, 0:phy_spo_lay%jmax-1, 0:phy_spo_lay%kmax-1)
                              ! $ u $ .     東西風速. Zonal wind
    real(DP), intent(in):: xyz_V (0:phy_spo_lay%imax-1, 0:phy_spo_lay%jmax-1, 0:phy_spo_lay%kmax-1)
                              ! $ v $ .     南北風速. Meridional wind
    real(DP), intent(in):: xy_Ps  (0:phy_spo_lay%imax-1, 0:phy_spo_lay%jmax-1)
                              ! $ p_s $ .   地表面気圧. Surface pressure

    real(DP), intent(out):: xyz_DUDtSpoDamping (0:phy_spo_lay%imax-1, 0:phy_spo_lay%jmax-1, 0:phy_spo_lay%kmax-1)
                              ! $ \DP{u}{t} $ . 
                              ! 減衰の効果による東西風速変化. 
                              ! Zonal wind tendency by damping effect
    real(DP), intent(out):: xyz_DVDtSpoDamping (0:phy_spo_lay%imax-1, 0:phy_spo_lay%jmax-1, 0:phy_spo_lay%kmax-1)
                              ! $ \DP{v}{t} $ . 
                              ! 減衰の効果による南北風速変化. 
                              ! Meridional wind tendency by damping effect
    logical, intent(in), optional:: historyput_flag
                              ! データ出力のフラグ. 
                              ! SetTime によって時刻を明示的に
                              ! 指定した場合には, この引数に
                              ! .true. または .false. を指定する
                              ! ことでデータ出力のオンオフを
                              ! 明示的に指定する必要があります. 
                              ! デフォルトは .false. です. 
                              ! 
                              ! Data output flag. 
                              ! When time is specified by "SetTime", 
                              ! explicit specification of data output 
                              ! on/off by specifying ".true." or ".false." 
                              ! to this argument. 
                              ! Default value is ".false.". 
                              ! 
    logical, intent(out), optional:: err
                              ! 例外処理用フラグ. 
                              ! デフォルトでは, この手続き内でエラーが
                              ! 生じた場合, プログラムは強制終了します. 
                              ! 引数 *err* が与えられる場合, 
                              ! プログラムは強制終了せず, 代わりに
                              ! *err* に .true. が代入されます. 
                              !
                              ! Exception handling flag. 
                              ! By default, when error occur in 
                              ! this procedure, the program aborts. 
                              ! If this *err* argument is given, 
                              ! .true. is substituted to *err* and 
                              ! the program does not abort. 

    !-----------------------------------
    !  *phy_spo_lay* から取り出される設定値
    !  Setting values fetched from *phy_spo_lay*
    integer:: imax            ! 経度格子点数. 
                              ! Number of grid points in longitude
    integer:: jmax            ! 緯度格子点数. 
                              ! Number of grid points in latitude
    integer:: kmax            ! 鉛直層数. 
                              ! Number of vertical level

    real(DP):: z_Sigma (0:phy_spo_lay%kmax-1)
                              ! $ \sigma $ レベル (整数). 
                              ! Full $ \sigma $ level

    real(DP):: DampingTime
                              ! 減衰の時定数 (単位: 秒). 
                              ! Damping time constant (units: seconds)
    real(DP):: UpperPress
                              ! 減衰範囲の上端の気圧. 
                              ! Pressure in the top within the range of damping
    real(DP):: LowerPress
                              ! 減衰範囲の下端の気圧. 
                              ! Pressure at bottom within the range of damping

    !-----------------------------------
    !  ヒストリファイルへのデータ出力設定
    !  Configure the settings for history data output
    character(STRING):: name = ''
                              ! 変数名. Variable identifier
    real:: time
                              ! 時刻. Time
    type(GT_HISTORY), pointer:: gthist =>null()
                              ! gt4_history モジュール用構造体. 
                              ! Derived type for "gt4_history" module

    !-----------------------------------
    !  作業変数
    !  Work variables
    real(DP):: xyz_Press (0:phy_spo_lay%imax-1, 0:phy_spo_lay%jmax-1, 0:phy_spo_lay%kmax-1)
                              ! $ p $ .     気圧. Air pressure
    real(DP):: yz_UAvrLon (0:phy_spo_lay%jmax-1, 0:phy_spo_lay%kmax-1)
                              ! 東西風速の東西平均
                              ! Zonal mean zonal wind
    real(DP):: yz_VAvrLon (0:phy_spo_lay%jmax-1, 0:phy_spo_lay%kmax-1)
                              ! 南北風速の東西平均
                              ! Zonal mean meritional wind

    integer:: i, j, k         ! DO ループ用作業変数
                              ! Work variables for DO loop

    integer:: stat
    character(STRING):: cause_c
    character(*), parameter:: subname = 'PhySpoLayDamping'
  continue
    call BeginSub( subname )
    stat = DC_NOERR
    cause_c = ''

    !-----------------------------------------------------------------
    !  初期設定のチェック
    !  Check initialization
    !-----------------------------------------------------------------
    if ( .not. phy_spo_lay % initialized ) then
      stat = DC_ENOTINIT
      cause_c = 'PHYSPOLAY'
      goto 999
    end if

    !-----------------------------------------------------------------
    !  *phy_spo_lay* に格納されている設定値の取り出し
    !  Fetch setting values stored in *phy_spo_lay*
    !-----------------------------------------------------------------
    imax  = phy_spo_lay % imax 
    jmax  = phy_spo_lay % jmax 
    kmax  = phy_spo_lay % kmax 

    z_Sigma  =  phy_spo_lay % z_Sigma 

    DampingTime = phy_spo_lay % DampingTime 
    UpperPress  = phy_spo_lay % UpperPress 
    LowerPress  = phy_spo_lay % LowerPress 

    !-----------------------------------------------------------------
    !  時間変化の演算
    !  Calculate tendency
    !-----------------------------------------------------------------
    do k = 0, kmax-1
      xyz_Press(:,:,k) = xy_Ps * z_Sigma(k)
    end do

    yz_UAvrLon = ya_AvrLon_xya( xyz_U, phy_spo_lay % intavr_opr )
    yz_VAvrLon = ya_AvrLon_xya( xyz_V, phy_spo_lay % intavr_opr )

    xyz_DUDtSpoDamping = 0.0_DP
    xyz_DVDtSpoDamping = 0.0_DP

    do k = 0, kmax-1
      do j = 0, jmax-1
        do i = 0, imax-1
          if ( ( xyz_Press(i,j,k) < LowerPress ) .and. ( xyz_Press(i,j,k) > UpperPress ) ) then
            xyz_DUDtSpoDamping(i,j,k) = ( yz_UAvrLon(j,k) - xyz_U(i,j,k) ) / DampingTime
            xyz_DVDtSpoDamping(i,j,k) = ( yz_VAvrLon(j,k) - xyz_V(i,j,k) ) / DampingTime
          end if
        end do
      end do
    end do

    !----------------------------------------------------------------
    !  ヒストリファイルへのデータ出力
    !  History data output
    !----------------------------------------------------------------

    !-----------------------------------
    !  xyz_DUDtSpoDamping の出力
    !  Output "xyz_DUDtSpoDamping"
    name = 'DUDtSpoDamping'

    ! 出力のチェック. 
    !   * gthist (gt4_history#GT_HISTORY), time (単精度実数型) が設定される.
    ! Check for output.
    !   * "gthist" (gt4_history#GT_HISTORY), time (real) are configured.
    call output_check ! これは内部サブルーチン. This is an internal subroutine

    ! 出力データを引数 array に渡す.
    ! Give output data to argument "array"
    if ( associated( gthist ) ) then
      call HistoryPut( history = gthist, varname = name, array = xyz_DUDtSpoDamping, time = time, quiet = .false., err = err )                        ! (out)
    end if

    !-----------------------------------
    !  xyz_DVDtSpoDamping の出力
    !  Output "xyz_DVDtSpoDamping"
    name = 'DVDtSpoDamping'

    ! 出力のチェック. 
    !   * gthist (gt4_history#GT_HISTORY), time (単精度実数型) が設定される.
    ! Check for output.
    !   * "gthist" (gt4_history#GT_HISTORY), time (real) are configured.
    call output_check ! これは内部サブルーチン. This is an internal subroutine

    ! 出力データを引数 array に渡す.
    ! Give output data to argument "array"
    if ( associated( gthist ) ) then
      call HistoryPut( history = gthist, varname = name, array = xyz_DVDtSpoDamping, time = time, quiet = .false., err = err )                        ! (out)
    end if

!!$    !-----------------------------------
!!$    !  y_Data2 の出力
!!$    !  Output "y_Data2"
!!$    name = 'Data2'
!!$
!!$    ! 出力のチェック. 
!!$    !   * gthist (gt4_history#GT_HISTORY), time (単精度実数型) が設定される.
!!$    ! Check for output.
!!$    !   * "gthist" (gt4_history#GT_HISTORY), time (real) are configured.
!!$    call output_check ! これは内部サブルーチン. This is an internal subroutine
!!$
!!$    ! 出力データを引数 array に渡す.
!!$    ! Give output data to argument "array"
!!$    if ( associated( gthist ) ) then
!!$      call HistoryPut( &
!!$        & history = gthist, &                ! (inout)
!!$        & varname = name, array = y_Data2, & ! (in)
!!$        & time = time, quiet = .false., &    ! (in)
!!$        & err = err )                        ! (out)
!!$    end if

    !-----------------------------------------------------------------
    !  時刻の更新
    !  Update time
    !-----------------------------------------------------------------
    phy_spo_lay % current_time = phy_spo_lay % current_time + phy_spo_lay % delta_time

    !-----------------------------------------------------------------
    !  終了処理, 例外処理
    !  Termination and Exception handling
    !-----------------------------------------------------------------
999 continue
    call StoreError( stat, subname, err, cause_c )
    call EndSub( subname )

  contains

    subroutine output_check
      !
      ! 変数 *name* を出力するかどうかをチェックします. 
      ! 出力に関する情報は dcm_sam_cod % gthstnml から取り出されます. 
      !
      ! 変数 *name* に関して出力するよう設定されている場合には, 
      ! *gthist* に出力先ファイルの gt4_history#GT_HISTORY
      ! 型変数を結合させます. そうでない場合は, *gthist* を空状態にします. 
      !
      ! また, 現在時刻を *time* に設定します. 
      !
      ! Check whether to output variable *name*. 
      ! Information about output is taken out of "dcm_sam_cod % gthstnml". 
      !
      ! When output is done for the variable *name*, *gthist* is 
      ! associated with "gt4_history#GT_HISTORY" variable of
      ! the output file. Otherwise, *gthist* is nullified. 
      !
      ! Moreover, current time is set to *time*. 
      !
      character(TOKEN):: interval_unit
                                ! ヒストリデータの出力間隔の単位. 
                                ! Unit for interval of history data output
    continue

      nullify( gthist )
      time = 0.0

      if ( HstNmlInfoOutputValid( phy_spo_lay % gthstnml, name ) ) then
        call HstNmlInfoInquire( gthstnml = phy_spo_lay % gthstnml, name = name, interval_unit = interval_unit )          ! (out)

        time = real( EvalbyUnit( phy_spo_lay % current_time, interval_unit ) )
        if ( present_and_true( historyput_flag ) ) time = 0.0

        call HstNmlInfoAssocGtHist( gthstnml = phy_spo_lay % gthstnml, name = name, history = gthist, err = err )            ! (out)
        if ( present_and_true( err ) ) then
          nullify( gthist )
          return
        end if

        if ( .not. HistoryInitialized( gthist ) ) nullify( gthist )
      end if

    end subroutine output_check
  end subroutine PhySpoLayDamping
Subroutine :
nmlfile :character(*), intent(in)
: NAMELIST ファイルの名称. NAMELIST file name
DampingTime :real(DP), intent(inout)
: 減衰の時定数 (単位: 秒). Damping time constant (units: seconds)
UpperPress :real(DP), intent(inout)
: 減衰範囲の上端の気圧. Pressure in the top within the range of damping
LowerPress :real(DP), intent(inout)
: 減衰範囲の下端の気圧. Pressure at bottom within the range of damping
gthstnml :type(GTHST_NMLINFO), intent(inout)
: NAMELIST#phy_sponge_layer_history_nml から入手される個別のデータ出力情報.

初期設定やデフォルト値の設定などを 行った後に与えること.

Individual data output information from "NAMELIST#phy_sponge_layer_history_nml".

Before this argument is given to this procedure, initialize and configure the defaut settings.

err :logical, intent(out), optional
: 例外処理用フラグ. デフォルトでは, この手続き内でエラーが 生じた場合, プログラムは強制終了します. 引数 err が与えられる場合, プログラムは強制終了せず, 代わりに err に .true. が代入されます.

Exception handling flag. By default, when error occur in this procedure, the program aborts. If this err argument is given, .true. is substituted to err and the program does not abort.

NAMELIST ファイル nmlfile から値を入力するための サブルーチンです. PhySpoLayCreate 内で呼び出されることを想定しています.

値が NAMELIST ファイル内で指定されていない場合には, 入力された値がそのまま返ります.

なお, nmlfile に空文字が与えられた場合, または 与えられた nmlfile を読み込むことができない場合, プログラムはエラーを発生させます.

This is a subroutine to input values from NAMELIST file nmlfile. This subroutine is expected to be called by "PhySpoLayCreate".

A value not specified in NAMELIST file is returned without change.

If nmlfile is empty, or nmlfile can not be read, error is occurred.

This procedure input/output NAMELIST#phy_sponge_layer_nml, NAMELIST#phy_sponge_layer_history_nml .

[Source]

  subroutine PhySpoLayNmlRead( nmlfile, DampingTime, UpperPress, LowerPress, gthstnml, err )
    !
    ! NAMELIST ファイル *nmlfile* から値を入力するための
    ! サブルーチンです. PhySpoLayCreate 
    ! 内で呼び出されることを想定しています. 
    !
    ! 値が NAMELIST ファイル内で指定されていない場合には, 
    ! 入力された値がそのまま返ります. 
    !
    ! なお, *nmlfile* に空文字が与えられた場合, または
    ! 与えられた *nmlfile* を読み込むことができない場合, 
    ! プログラムはエラーを発生させます. 
    !
    ! This is a subroutine to input values from 
    ! NAMELIST file *nmlfile*. This subroutine is expected to be 
    ! called by "PhySpoLayCreate". 
    !
    ! A value not specified in NAMELIST file is returned 
    ! without change. 
    !
    ! If *nmlfile* is empty, or *nmlfile* can not be read, 
    ! error is occurred. 
    !
    use dc_trace, only: BeginSub, EndSub
    use dc_string, only: PutLine, Printf, Split, StrInclude, StoA, JoinChar
    use dc_types, only: DP, STRING, TOKEN, STDOUT
    use dc_iounit, only: FileOpen
    use dc_message, only: MessageNotify
    use dc_present, only: present_and_true
    use dc_date, only: DCDiffTimeCreate
    use dc_error, only: StoreError, DC_NOERR, DC_ENOFILEREAD, DC_ENOTINIT
    use gt4_history_nmlinfo, only: GTHST_NMLINFO, HstNmlInfoAdd, HstNmlInfoInquire, HstNmlInfoInitialized, HstNmlInfoPutLine
    implicit none
    character(*), intent(in):: nmlfile
                              ! NAMELIST ファイルの名称. 
                              ! NAMELIST file name

    real(DP), intent(inout):: DampingTime
                              ! 減衰の時定数 (単位: 秒). 
                              ! Damping time constant (units: seconds)
    real(DP), intent(inout):: UpperPress
                              ! 減衰範囲の上端の気圧. 
                              ! Pressure in the top within the range of damping
    real(DP), intent(inout):: LowerPress
                              ! 減衰範囲の下端の気圧. 
                              ! Pressure at bottom within the range of damping
    type(GTHST_NMLINFO), intent(inout):: gthstnml
                              ! NAMELIST#phy_sponge_layer_history_nml 
                              ! から入手される個別のデータ出力情報. 
                              ! 
                              ! 初期設定やデフォルト値の設定などを
                              ! 行った後に与えること. 
                              ! 
                              ! Individual data output information from 
                              ! "NAMELIST#phy_sponge_layer_history_nml". 
                              ! 
                              ! Before this argument is given to 
                              ! this procedure, initialize and 
                              ! configure the defaut settings. 
    logical, intent(out), optional:: err
                              ! 例外処理用フラグ. 
                              ! デフォルトでは, この手続き内でエラーが
                              ! 生じた場合, プログラムは強制終了します. 
                              ! 引数 *err* が与えられる場合, 
                              ! プログラムは強制終了せず, 代わりに
                              ! *err* に .true. が代入されます. 
                              !
                              ! Exception handling flag. 
                              ! By default, when error occur in 
                              ! this procedure, the program aborts. 
                              ! If this *err* argument is given, 
                              ! .true. is substituted to *err* and 
                              ! the program does not abort. 

    namelist /phy_sponge_layer_nml/ DampingTime, UpperPress, LowerPress
                              ! phy_sponge_layer モジュール用
                              ! NAMELIST 変数群名. 
                              !
                              ! phy_sponge_layer#PhySpoLayCreate 
                              ! を使用する際に, オプショナル引数 *nmlfile* 
                              ! へ NAMELIST ファイル名を指定することで, 
                              ! そのファイルからこの NAMELIST 変数群を
                              ! 読み込みます. 
                              !
                              ! NAMELIST group name for 
                              ! "phy_sponge_layer" module. 
                              ! 
                              ! If a NAMELIST filename is specified to 
                              ! an optional argument *nmlfile* when 
                              ! "phy_sponge_layer#PhySpoLayCreate" 
                              ! is used, this NAMELIST group is 
                              ! loaded from the file. 

    character(STRING):: name
                              ! 変数名. 
                              ! 空白の場合には, この他の設定値は
                              ! phy_sponge_layer モジュールにおいて
                              ! 出力されるデータ全ての
                              ! デフォルト値となります. 
                              ! 
                              ! "Data1,Data2" のようにカンマで区切って複数
                              ! の変数を指定することも可能です. 
                              ! 
                              ! Variable identifier. 
                              ! If blank is given, other values are 
                              ! used as default values of output data 
                              ! in "phy_sponge_layer". 
                              ! 
                              ! Multiple variables can be specified 
                              ! as "Data1,Data2" too. Delimiter is comma. 
    character(STRING):: file
                              ! 出力ファイル名. 
                              ! これはデフォルト値としては使用されません. 
                              ! *name* に値が設定されている時のみ有効です. 
                              ! 
                              ! Output file name. 
                              ! This is not used as default value. 
                              ! This value is valid only when *name* is 
                              ! specified. 

    real:: interval_value
                              ! ヒストリデータの出力間隔の数値. 
                              ! 負の値を与えると, 出力を抑止します. 
                              ! Numerical value for interval of history data output
                              ! Negative values suppresses output.
    character(TOKEN):: interval_unit
                              ! ヒストリデータの出力間隔の単位. 
                              ! Unit for interval of history data output
    character(TOKEN):: precision
                              ! ヒストリデータの精度. 
                              ! Precision of history data
    logical:: average
                              ! 出力データの平均化フラグ. 
                              ! Flag for average of output data
    character(STRING):: fileprefix
                              ! ヒストリデータのファイル名の接頭詞. 
                              ! Prefixes of history data filenames

    namelist /phy_sponge_layer_history_nml/ name, file, interval_value, interval_unit, precision, fileprefix, average
                              ! phy_sponge_layer モジュールのヒストリデータ用
                              ! NAMELIST 変数群名. 
                              !
                              ! phy_sponge_layer#PhySpoLayCreate 
                              ! を使用する際に, オプショナル引数 *nmlfile* 
                              ! へ NAMELIST ファイル名を指定することで, 
                              ! そのファイルからこの NAMELIST 変数群を
                              ! 読み込みます. 
                              !
                              ! NAMELIST group name for 
                              ! history data of "phy_sponge_layer" module. 
                              ! 
                              ! If a NAMELIST filename is specified to 
                              ! an optional argument *nmlfile* when 
                              ! "phy_sponge_layer#PhySpoLayCreate" 
                              ! is used, this NAMELIST group is 
                              ! loaded from the file. 

    !-----------------------------------
    !  作業変数
    !  Work variables
    integer:: stat
    character(STRING):: cause_c
    integer:: unit_nml        ! NAMELIST ファイルオープン用装置番号. 
                              ! Unit number for NAMELIST file open
    integer:: iostat_nml      ! NAMELIST 読み込み時の IOSTAT. 
                              ! IOSTAT of NAMELIST read
    character(TOKEN):: pos_nml
                              ! NAMELIST 読み込み時のファイル位置. 
                              ! File position of NAMELIST read
    character(*), parameter:: subname = 'PhySpoLayNmlRead'
  continue
    call BeginSub( subname )
    stat = DC_NOERR
    cause_c = ''

    !-----------------------------------------------------------------
    !  初期設定のチェック
    !  Check initialization
    !-----------------------------------------------------------------

    !-----------------------------------------------------------------
    !  文字型引数を NAMELIST 変数群へ代入
    !  Substitute character arguments to NAMELIST group
    !-----------------------------------------------------------------
!!$    key00       = key00_

    !----------------------------------------------------------------
    !  NAMELIST ファイルのオープン
    !  Open NAMELIST file
    !----------------------------------------------------------------
    call FileOpen( unit = unit_nml, file = nmlfile, mode = 'r', err = err )                   ! (out)
    if ( present_and_true(err) ) then
      stat = DC_ENOFILEREAD
      cause_c = nmlfile
      goto 999
    end if

    !-----------------------------------------------------------------
    !  NAMELIST 変数群の取得
    !  Get NAMELIST group
    !-----------------------------------------------------------------

    !-------------------------
    !  係数などの取得
    !  Get coefficients etc.
    rewind( unit_nml )
    read( unit = unit_nml, nml = phy_sponge_layer_nml, iostat = iostat_nml )        ! (out)
    if ( iostat_nml == 0 ) then
      call MessageNotify( 'M', subname, 'NAMELIST group "%c" is loaded from "%c".', c1 = 'phy_sponge_layer_nml', c2 = trim(nmlfile) )
      write(STDOUT, nml = phy_sponge_layer_nml)
    else
      call MessageNotify( 'W', subname, 'NAMELIST group "%c" is not found in "%c" (iostat=%d).', c1 = 'phy_sponge_layer_nml', c2 = trim(nmlfile), i = (/iostat_nml/) )
    end if

    !-------------------------
    !  出力データの個別情報の取得
    !  Get individual information of output data
    rewind( unit_nml )
    iostat_nml = 0
    pos_nml = ''
    do while ( trim(pos_nml) /= 'APPEND' .and. iostat_nml == 0 )

      name                   = ''
      file                   = ''
      call HstNmlInfoInquire( gthstnml = gthstnml, interval_value = interval_value, interval_unit = interval_unit, precision = precision, average = average, fileprefix = fileprefix )          ! (out)

      read( unit = unit_nml, nml = phy_sponge_layer_history_nml, iostat = iostat_nml )                ! (out)
      inquire( unit = unit_nml, position = pos_nml )   ! (out)

      if ( iostat_nml == 0 ) then
        call MessageNotify( 'M', subname, 'NAMELIST group "%c" is loaded from "%c".', c1='phy_sponge_layer_history_nml', c2=trim(nmlfile) )
        write(STDOUT, nml = phy_sponge_layer_history_nml)

        call HstNmlInfoAdd( gthstnml = gthstnml, name = name, file = file, interval_value = interval_value, interval_unit = interval_unit, precision = precision, average = average, fileprefix = fileprefix )          ! (in)

      else
        call MessageNotify( 'W', subname, 'NAMELIST group "%c" is not found in "%c" any more (iostat=%d).', c1='phy_sponge_layer_history_nml', c2=trim(nmlfile), i = (/iostat_nml/) )
      end if
    end do

    close( unit_nml )

    !-----------------------------------------------------------------
    !  NAMELIST 変数群を文字型引数へ代入
    !  Substitute NAMELIST group to character arguments
    !-----------------------------------------------------------------
!!$    key00_       = key00

    !-----------------------------------------------------------------
    !  終了処理, 例外処理
    !  Termination and Exception handling
    !-----------------------------------------------------------------
999 continue
    call StoreError( stat, subname, err, cause_c )
    call EndSub( subname )

  end subroutine PhySpoLayNmlRead
Subroutine :

変数 name を出力するかどうかをチェックします. 出力に関する情報は dcm_sam_cod % gthstnml から取り出されます.

変数 name に関して出力するよう設定されている場合には, gthist に出力先ファイルの gt4_history#GT_HISTORY 型変数を結合させます. そうでない場合は, gthist を空状態にします.

また, 現在時刻を time に設定します.

Check whether to output variable name. Information about output is taken out of "dcm_sam_cod % gthstnml".

When output is done for the variable name, gthist is associated with "gt4_history#GT_HISTORY" variable of the output file. Otherwise, gthist is nullified.

Moreover, current time is set to time.

[Source]

    subroutine output_check
      !
      ! 変数 *name* を出力するかどうかをチェックします. 
      ! 出力に関する情報は dcm_sam_cod % gthstnml から取り出されます. 
      !
      ! 変数 *name* に関して出力するよう設定されている場合には, 
      ! *gthist* に出力先ファイルの gt4_history#GT_HISTORY
      ! 型変数を結合させます. そうでない場合は, *gthist* を空状態にします. 
      !
      ! また, 現在時刻を *time* に設定します. 
      !
      ! Check whether to output variable *name*. 
      ! Information about output is taken out of "dcm_sam_cod % gthstnml". 
      !
      ! When output is done for the variable *name*, *gthist* is 
      ! associated with "gt4_history#GT_HISTORY" variable of
      ! the output file. Otherwise, *gthist* is nullified. 
      !
      ! Moreover, current time is set to *time*. 
      !
      character(TOKEN):: interval_unit
                                ! ヒストリデータの出力間隔の単位. 
                                ! Unit for interval of history data output
    continue

      nullify( gthist )
      time = 0.0

      if ( HstNmlInfoOutputValid( phy_spo_lay % gthstnml, name ) ) then
        call HstNmlInfoInquire( gthstnml = phy_spo_lay % gthstnml, name = name, interval_unit = interval_unit )          ! (out)

        time = real( EvalbyUnit( phy_spo_lay % current_time, interval_unit ) )
        if ( present_and_true( historyput_flag ) ) time = 0.0

        call HstNmlInfoAssocGtHist( gthstnml = phy_spo_lay % gthstnml, name = name, history = gthist, err = err )            ! (out)
        if ( present_and_true( err ) ) then
          nullify( gthist )
          return
        end if

        if ( .not. HistoryInitialized( gthist ) ) nullify( gthist )
      end if

    end subroutine output_check
Subroutine :

変数 name に関して出力ファイルの初期設定を行います. 出力ファイル名や出力間隔などの情報は phy_spo_lay % gthstnml から取り出されます.

変数 name に関して出力が行われる場合には, gthist に出力先ファイルの gt4_history#GT_HISTORY 型変数を結合させます. そうでない場合は, gthist を空状態にします.

また, 出力データの精度を precision に, 出力データ平均化の可否を average に設定します.

標準出力に表示される変数リスト registered_varnamesname, longname, dims, units が登録されます.

An output file is initialized for a variable name. Information such as the output filename and output intervals is taken out of "phy_spo_lay % gthstnml".

When output is done for the variable name, gthist is associated with the "gt4_history#GT_HISTORY" variable of the output file. Otherwise, gthist is nullified.

Moreover, the accuracy of output data is set to precision, and right or wrong of averaging the output data is set to average.

name, longname, dims, units are registered to a list of variables registered_varnames that is printed to standard output.

[Source]

    subroutine output_init
      !
      ! 変数 *name* に関して出力ファイルの初期設定を行います. 
      ! 出力ファイル名や出力間隔などの情報は phy_spo_lay % gthstnml 
      ! から取り出されます. 
      ! 
      ! 変数 *name* に関して出力が行われる場合には, 
      ! *gthist* に出力先ファイルの gt4_history#GT_HISTORY
      ! 型変数を結合させます. そうでない場合は, *gthist* を空状態にします. 
      !
      ! また, 出力データの精度を precision に, 
      ! 出力データ平均化の可否を average に設定します. 
      ! 
      ! 標準出力に表示される変数リスト *registered_varnames* に
      ! *name*, *longname*, *dims*, *units* が登録されます. 
      !
      ! An output file is initialized for a variable *name*. 
      ! Information such as the output filename and output intervals 
      ! is taken out of "phy_spo_lay % gthstnml". 
      !
      ! When output is done for the variable *name*, *gthist* is 
      ! associated with the "gt4_history#GT_HISTORY" variable of 
      ! the output file. Otherwise, *gthist* is nullified. 
      !
      ! Moreover, the accuracy of output data is set to *precision*, and 
      ! right or wrong of averaging the output data is set to *average*. 
      !
      ! *name*, *longname*, *dims*, *units* are registered to 
      ! a list of variables *registered_varnames* that is printed to 
      ! standard output. 
      !
      use dc_date, only: DCDiffTimeCreate, EvalSec, EvalByUnit
      use gt4_history_nmlinfo, only: GTHST_NMLINFO, HstNmlInfoInitialized, HstNmlInfoInquire, HstNmlInfoOutputValid, HstNmlInfoAssocGtHist, HstNmlInfoPutLine, HstNmlInfoSetValidName
      use gt4_history, only: GT_HISTORY, HistoryCreate, HistoryAddVariable, HistoryPut, HistoryAddAttr, HistoryInitialized

      !-----------------------------------
      !  作業変数
      !  Work variables
      character(STRING):: file
                                ! ヒストリデータのファイル名. 
                                ! History data filenames
      character(STRING):: dims_str
                                ! 座標軸のリスト. 
                                ! List of axes
      real:: interval_value
                                ! ヒストリデータの出力間隔の数値. 
                                ! Numerical value for interval of history data output
      character(TOKEN):: interval_unit
                                ! ヒストリデータの出力間隔の単位. 
                                ! Unit for interval of history data output
    continue
      !-----------------------------------------------------------------
      !  標準出力に表示される変数の登録
      !  Register a variable name for print to standard output
      !-----------------------------------------------------------------
      if ( allocated(dims) ) then
        dims_str = JoinChar( dims, ',' )
      else
        dims_str = ''
      end if
      call DCHashPut( hashv = registered_varnames, key = name, value = trim( longname ) // ' [' // trim( units ) // '] {' // trim( dims_str ) // '}' ) ! (in)

      !-----------------------------------------------------------------
      !  変数の初期化
      !  Initialize variable
      !-----------------------------------------------------------------
      nullify( gthist )
      precision = 'float'
      average = .false.

      !-----------------------------------------------------------------
      !  変数名の有効性を設定
      !  Set validation of the variable name
      !-----------------------------------------------------------------
      call HstNmlInfoSetValidName( gthstnml = phy_spo_lay % gthstnml, name = name )                        ! (in)

      !-----------------------------------------------------------------
      !  出力が有効かどうかを確認する
      !  Confirm whether the output is effective
      !-----------------------------------------------------------------
      if ( .not. HstNmlInfoOutputValid( phy_spo_lay % gthstnml, name ) ) then
        return
      end if

      !-----------------------------------------------------------------
      !  GT_HISTORY 変数の取得
      !  Get "GT_HISTORY" variable
      !-----------------------------------------------------------------
      call HstNmlInfoAssocGtHist( gthstnml = phy_spo_lay % gthstnml, name = name, history = gthist, err = err )                          ! (out)
      if ( present_and_true( err ) ) return

      call HstNmlInfoInquire( gthstnml = phy_spo_lay % gthstnml, name = name, precision = precision, average = average, err = err )                          ! (out)
      if ( present_and_true( err ) ) return

      !-----------------------------------------------------------------
      !  GT_HISTORY 変数の初期設定の確認
      !  Check initialization of "GT_HISTORY" variable
      !-----------------------------------------------------------------
      if ( HistoryInitialized( gthist ) ) then

        !---------------------------------------------------------------
        !  HistoryAddVariable による変数作成
        !  A variable is created by "HistoryAddVariable"
        !---------------------------------------------------------------
        call HistoryAddVariable( history = gthist, varname = name,         dims = dims, longname = longname,    units = units, xtype = precision, average = average )  ! (in)
        return
      end if

      !-----------------------------------------------------------------
      !  HistoryCreate のための設定値の取得
      !  Get the settings for "HistoryCreate"
      !-----------------------------------------------------------------
      call HstNmlInfoInquire( gthstnml = phy_spo_lay % gthstnml, name = name, file = file, interval_unit = interval_unit, interval_value = interval_value, err = err )                          ! (out)
      if ( present_and_true( err ) ) return

      !-----------------------------------------------------------------
      !  HistoryCreate によるファイル作成
      !  Files are created by "HistoryCreate"
      !-----------------------------------------------------------------
      call HistoryCreate( history = gthist, file = file, title = 'Sponge layer for damping velocity', source = 'dcpam4 project : ' // trim(version), institution = 'GFD Dennou Club', dims = StoA( 'lon', 'lat', 'sig', 'time' ), dimsizes = (/ phy_spo_lay % imax, phy_spo_lay % jmax, phy_spo_lay % kmax, 0 /), longnames = StoA('longitude', 'latitude', 'sigma at layer midpoints', 'time'), units = StoA( 'degree_east', 'degree_north', '1', interval_unit ), origin = real( EvalbyUnit( phy_spo_lay % current_time, interval_unit) ), interval = interval_value, err = err )                                            ! (out)
      if ( present_and_true( err ) ) then
        nullify( gthist )
        return
      end if

      call HistoryAddAttr( history = gthist, varname = 'lon', attrname = 'standard_name', value = 'longitude' )                          ! (in)
      call HistoryAddAttr( history = gthist, varname = 'lat', attrname = 'standard_name', value = 'latitude' )                           ! (in)
      call HistoryAddAttr( history = gthist, varname = 'sig', attrname = 'standard_name', value = 'atmosphere_sigma_coordinate' )          ! (in)
      call HistoryAddAttr( history = gthist, varname = 'time', attrname = 'standard_name', value = 'time' )                                ! (in)
      call HistoryAddAttr( history = gthist, varname = 'sig', attrname = 'positive', value = 'down' )                                 ! (in)

      call HistoryPut( history = gthist, varname = 'lon', array = phy_spo_lay % x_Lon / PI * 180.0_DP ) ! (in)
      call HistoryPut( history = gthist, varname = 'lat', array = phy_spo_lay % y_Lat / PI * 180.0_DP ) ! (in)
      call HistoryPut( history = gthist, varname = 'sig', array = phy_spo_lay % z_Sigma ) ! (in)

      call HistoryAddVariable( history = gthist, varname = 'lon_weight', dims = StoA('lon'), longname = 'weight for integration in longitude', units = 'radian', xtype = 'double' )                ! (in)
      call HistoryAddAttr( history = gthist, varname = 'lon', attrname = 'gt_calc_weight', value = 'lon_weight' )                              ! (in)
      call HistoryPut( history = gthist, varname = 'lon_weight', array = phy_spo_lay % x_Lon_Weight )    ! (in)

      call HistoryAddVariable( history = gthist, varname = 'lat_weight', dims = StoA('lat'), longname = 'weight for integration in latitude', units = 'radian', xtype = 'double' )                ! (in)
      call HistoryAddAttr( history = gthist, varname = 'lat', attrname = 'gt_calc_weight', value = 'lat_weight' )                              ! (in)
      call HistoryPut( history = gthist, varname = 'lat_weight', array = phy_spo_lay % y_Lat_Weight )    ! (in)

      call HistoryAddVariable( history = gthist, varname = 'sig_weight', dims = StoA('sig'), longname = 'weight for integration in sigma', units = '1', xtype = 'float' )                      ! (in)
      call HistoryAddAttr( history = gthist, varname = 'sig', attrname = 'gt_calc_weight', value = 'sig_weight' )                              ! (in)
      call HistoryPut( history = gthist, varname = 'sig_weight', array = phy_spo_lay % z_DelSigma )      ! (in)

      !-----------------------------------------------------------------
      !  HistoryAddVariable による変数作成
      !  A variable is created by "HistoryAddVariable"
      !-----------------------------------------------------------------
      if ( HistoryInitialized( gthist ) ) then
        call HistoryAddVariable( history = gthist, varname = name,         dims = dims, longname = longname,    units = units, xtype = precision, average = average )  ! (in)
      else
        nullify( gthist )
      end if

    end subroutine output_init
version
Constant :
version = ’$Name: dcpam4-20080609-1 $’ // ’$Id: phy_sponge_layer.f90,v 1.3 2008-06-01 16:06:53 morikawa Exp $’ :character(*), parameter

[Validate]