Class timeset
In: setup/timeset.f90

時刻管理

Time control

Note that Japanese and English are described in parallel.

時刻情報の管理を行います.

Information of time is controlled.

Variables List

Nstep :総ステップ数
Cstep :現在のステップ数
DelTime :$ Delta t $ [s]
———— :————
Nstep :Number of full steps
Cstep :Current steps
DelTime :$ Delta t $ [s]

Procedures List

TimesetInit :timeset モジュールの初期化
TimesetClose :timeset モジュールの終了処理
TimesetProgress :時刻の進行
TimesetClockStart :計算時間計測開始
TimesetClockStop :計算時間計一時停止
TimesetIntStepEval :ステップ間隔見積り
TimesetGetStartTime :計算開始時刻の取得
TimesetGetEndTime :計算終了時刻の取得
TimesetGetDelTime :$ Delta t $ の取得
TimesetGetCurrentTime :現在時刻の取得
———— :————
TimesetInit :Initialize "timeset" module
TimesetClose :Terminate "timeset" module
TimesetProgress :Progress time
TimesetClockStart :Start measurement of computation time
TimesetClockStop :Pause measurement of computation time
TimesetIntStepEval :Eval interval of step
TimesetGetStartTime :Get start time of calculation
TimesetGetEndTime :Get end time of calculation
TimesetGetDelTime :Get $ Delta t $
TimesetGetCurrentTime :Get current time

NAMELIST

NAMELIST#timeset_nml

Methods

Included Modules

dc_types dc_message dc_date_types dc_clock namelist_util dc_iounit dc_date

Public Instance methods

Cstep
Variable :
Cstep :integer, save, public
: 現在のステップ数. Current steps
DelTime
Variable :
DelTime :real(DP), save, public
: $ Delta t $ [s]
Nstep
Variable :
Nstep :integer, save, public
: 総ステップ数. Number of full steps
Subroutine :
name :character(*), intent(in)
: モジュールの名称. Name of module

プログラム単位 (主にモジュールを想定) ごとの時間計測を開始します.

Start measurement of computation time by program unit (expected modules).

[Source]

  subroutine TimesetClockStart( name )
    !
    ! プログラム単位 (主にモジュールを想定) ごとの時間計測を開始します. 
    !
    ! Start measurement of computation time by program unit
    ! (expected modules). 

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

    ! CPU 時間計測
    ! CPU time monitor
    !
    use dc_clock, only: DCClockCreate, DCClockStart

    ! 宣言文 ; Declaration statements
    !
    implicit none
    character(*), intent(in):: name
                              ! モジュールの名称. 
                              ! Name of module

    ! 作業変数
    ! Work variables
    !
    integer:: i               ! clocks, clocks_name 用 DO ループ用作業変数
                              ! Work variables for DO loop for "clocks", "clocks_name"

    ! 実行文 ; Executable statement
    !

    if ( .not. CpuTimeMoniter ) return

    do i = 1, clk_proc_num
      if ( trim(clocks_name(i)) == trim(name) ) then
        call DCClockStart( clocks(i) ) ! (in)
        return
      end if
    end do

    call DCClockCreate( clocks(clk_proc_num + 1), name ) ! (in)
    call DCClockStart( clocks(clk_proc_num + 1) ) ! (in)
    clocks_name(clk_proc_num + 1) = name
    clk_proc_num = clk_proc_num + 1

  end subroutine TimesetClockStart
Subroutine :
name :character(*), intent(in)
: モジュールの名称. Name of module

プログラム単位 (主にモジュールを想定) ごとの時間計測を一時停止します.

Pause measurement of computation time by program unit (expected modules).

[Source]

  subroutine TimesetClockStop( name )
    !
    ! プログラム単位 (主にモジュールを想定) ごとの時間計測を一時停止します.
    !
    ! Pause measurement of computation time by program unit
    ! (expected modules). 

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

    ! CPU 時間計測
    ! CPU time monitor
    !
    use dc_clock, only: DCClockStop

    ! 宣言文 ; Declaration statements
    !
    implicit none
    character(*), intent(in):: name
                              ! モジュールの名称. 
                              ! Name of module

    ! 作業変数
    ! Work variables
    !
    integer:: i               ! clocks, clocks_name 用 DO ループ用作業変数
                              ! Work variables for DO loop for "clocks", "clocks_name"

    ! 実行文 ; Executable statement
    !

    if ( .not. CpuTimeMoniter ) return

    do i = 1, clk_proc_num
      if ( trim(clocks_name(i)) == trim(name) ) then
        call DCClockStop( clocks(i) ) ! (in)
        return
      end if
    end do

    call MessageNotify( 'W', module_name, ' name "%c" is not found in "TimesetClockStop"', c1 = trim(name) )

  end subroutine TimesetClockStop
Subroutine :

計算時間の総計を表示します.

Total computation time is printed.

[Source]

  subroutine TimesetClose
    !
    ! 計算時間の総計を表示します. 
    !
    ! Total computation time is printed. 

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

    ! CPU 時間計測
    ! CPU time monitor
    !
    use dc_clock, only: DCClockStop, DCClockResult, DCClockSetName, operator(+), operator(-)

    ! 宣言文 ; Declaration statements
    !
    implicit none
    integer:: i               ! clocks, clocks_name 用 DO ループ用作業変数
                              ! Work variables for DO loop for "clocks", "clocks_name"

    ! 作業変数
    ! Work variables
    !

    ! 実行文 ; Executable statement
    !

    if ( .not. CpuTimeMoniter ) return

    ! CPU 時間計測終了 (モデル全体)
    ! Stop CPU time monitoring (for entire model)
    !
    call DCClockStop( clocks(1) ) ! (in)

    ! 「その他」の CPU 時間を算出
    ! Calculate CPU time of "Others"
    !
    clocks(clk_proc_num + 1) = clocks(1)
    do i = 2, clk_proc_num
      clocks(clk_proc_num + 1) = clocks(clk_proc_num + 1) - clocks(i)
    end do
    call DCClockSetName( clocks(clk_proc_num + 1), 'others' )

    ! CPU 時間の総計を表示
    ! Print total CPU time
    !
    call DCClockResult( clocks(2:clk_proc_num + 1), total_auto = .true. ) ! (in)

  end subroutine TimesetClose
Subroutine :
time :real(DP), intent(out)
: 現在時刻. Current time
unit :character(*), intent(in), optional
: 現在時刻の単位. Unit of current time

現在時刻 time に返します. デフォルトでは, NAMELIST#timeset_nmlDelTimeUnit に指定されるものがその時刻の単位です. オプショナル引数 unit に単位を指定することでその単位に 応じた時刻が返ります.

Current time is returned to "time". By default, "DelTimeUnit"in "NAMELIST#timeset_nml" is units of the time. If a unit is specified to an optional argument "unit", a time according to the unit is returned.

[Source]

  subroutine TimesetGetCurrentTime( time, unit )
    !
    ! 現在時刻 time に返します. 
    ! デフォルトでは, NAMELIST#timeset_nml の
    ! DelTimeUnit に指定されるものがその時刻の単位です. 
    ! オプショナル引数 unit に単位を指定することでその単位に
    ! 応じた時刻が返ります. 
    !
    ! Current time is returned to "time". 
    ! By default, "DelTimeUnit"in "NAMELIST#timeset_nml" is units of 
    ! the time. 
    ! If a unit is specified to an optional argument "unit", 
    ! a time according to the unit is returned. 
    !

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

    ! 日付および時刻の取り扱い
    ! Date and time handler
    !
    use dc_date_types, only: DC_DIFFTIME
                              ! 日時の差を表現するデータ型. 
                              ! Data type for difference about date and time
    use dc_date, only: DCDiffTimeCreate, EvalSec, EvalByUnit, mod, toChar, operator(*), operator(==), operator(<), operator(>), operator(/), operator(+), operator(-)

    ! 宣言文 ; Declaration statements
    !
    implicit none
    real(DP), intent(out):: time
                              ! 現在時刻. Current time
    character(*), intent(in), optional:: unit
                              ! 現在時刻の単位. 
                              ! Unit of current time

    ! 作業変数
    ! Work variables
    !

    ! 実行文 ; Executable statement
    !
    if ( .not. present(unit) ) then
      time = EvalByUnit( DCdiffCurrentTime, DelTimeUnit )
    else
      time = EvalByUnit( DCdiffCurrentTime, unit )
    end if

  end subroutine TimesetGetCurrentTime
Subroutine :
time :real(DP), intent(out)
: $ Delta t $
unit :character(*), intent(in), optional
: $ Delta t $ の単位. Unit of $ Delta t $

$ Delta t $ を time に返します. デフォルトでは, NAMELIST#timeset_nmlDelTimeUnit に指定されるものがその時刻の単位です. オプショナル引数 unit に単位を指定することでその単位に 応じた時刻が返ります.

$ Delta t $ is returned to "time". By default, "DelTimeUnit"in "NAMELIST#timeset_nml" is units of the time. If a unit is specified to an optional argument "unit", a time according to the unit is returned.

[Source]

  subroutine TimesetGetDelTime( time, unit )
    !
    ! $ \Delta t $ を time に返します. 
    ! デフォルトでは, NAMELIST#timeset_nml の
    ! DelTimeUnit に指定されるものがその時刻の単位です. 
    ! オプショナル引数 unit に単位を指定することでその単位に
    ! 応じた時刻が返ります. 
    !
    ! $ \Delta t $  is returned to "time". 
    ! By default, "DelTimeUnit"in "NAMELIST#timeset_nml" is units of 
    ! the time. 
    ! If a unit is specified to an optional argument "unit", 
    ! a time according to the unit is returned. 
    !

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

    ! 日付および時刻の取り扱い
    ! Date and time handler
    !
    use dc_date_types, only: DC_DIFFTIME
                              ! 日時の差を表現するデータ型. 
                              ! Data type for difference about date and time
    use dc_date, only: DCDiffTimeCreate, EvalByUnit

    ! 宣言文 ; Declaration statements
    !
    implicit none
    real(DP), intent(out):: time
                              ! $ \Delta t $ 
    character(*), intent(in), optional:: unit
                              ! $ \Delta t $ の単位. 
                              ! Unit of $ \Delta t $ 

    ! 作業変数
    ! Work variables
    !
    type(DC_DIFFTIME):: dcdiff_time

    ! 実行文 ; Executable statement
    !

    if ( .not. present(unit) ) then
      time = DelTimeValue
    else
      call DCDiffTimeCreate( dcdiff_time, DelTimeValue, DelTimeUnit )        ! (in)
      time = EvalByUnit( dcdiff_time, unit )
    end if

  end subroutine TimesetGetDelTime
Subroutine :
time :real(DP), intent(out)
: 計算終了時刻. End time of calculation
unit :character(*), intent(in), optional
: 計算終了時刻の単位. Unit of end time of calculation

計算終了時刻を time に返します. デフォルトでは, NAMELIST#timeset_nmlEndTimeUnit に指定されるものがその時刻の単位です. オプショナル引数 unit に単位を指定することでその単位に 応じた時刻が返ります.

End time of calculation is returned to "time". By default, "EndTimeUnit"in "NAMELIST#timeset_nml" is units of the time. If a unit is specified to an optional argument "unit", a time according to the unit is returned.

[Source]

  subroutine TimesetGetEndTime( time, unit )
    !
    ! 計算終了時刻を time に返します. 
    ! デフォルトでは, NAMELIST#timeset_nml の
    ! EndTimeUnit に指定されるものがその時刻の単位です. 
    ! オプショナル引数 unit に単位を指定することでその単位に
    ! 応じた時刻が返ります. 
    !
    ! End time of calculation is returned to "time". 
    ! By default, "EndTimeUnit"in "NAMELIST#timeset_nml" is units of 
    ! the time. 
    ! If a unit is specified to an optional argument "unit", 
    ! a time according to the unit is returned. 
    !

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

    ! 日付および時刻の取り扱い
    ! Date and time handler
    !
    use dc_date_types, only: DC_DIFFTIME
                              ! 日時の差を表現するデータ型. 
                              ! Data type for difference about date and time
    use dc_date, only: DCDiffTimeCreate, EvalByUnit

    ! 宣言文 ; Declaration statements
    !
    implicit none
    real(DP), intent(out):: time
                              ! 計算終了時刻. 
                              ! End time of calculation
    character(*), intent(in), optional:: unit
                              ! 計算終了時刻の単位. 
                              ! Unit of end time of calculation

    ! 作業変数
    ! Work variables
    !
    type(DC_DIFFTIME):: dcdiff_time

    ! 実行文 ; Executable statement
    !

    if ( .not. present(unit) ) then
      time = EndTimeValue
    else
      call DCDiffTimeCreate( dcdiff_time, EndTimeValue, EndTimeUnit )        ! (in)
      time = EvalByUnit( dcdiff_time, unit )
    end if

  end subroutine TimesetGetEndTime
Subroutine :
time :real(DP), intent(out)
: 計算開始時刻. Start time of calculation
unit :character(*), intent(in), optional
: 計算開始時刻の単位. Unit of start time of calculation

計算開始時刻を time に返します. デフォルトでは, NAMELIST#timeset_nmlStartTimeUnit に指定されるものがその時刻の単位です. オプショナル引数 unit に単位を指定することでその単位に 応じた時刻が返ります.

Start time of calculation is returned to "time". By default, "StartTimeUnit"in "NAMELIST#timeset_nml" is units of the time. If a unit is specified to an optional argument "unit", a time according to the unit is returned.

[Source]

  subroutine TimesetGetStartTime( time, unit )
    !
    ! 計算開始時刻を time に返します. 
    ! デフォルトでは, NAMELIST#timeset_nml の
    ! StartTimeUnit に指定されるものがその時刻の単位です. 
    ! オプショナル引数 unit に単位を指定することでその単位に
    ! 応じた時刻が返ります. 
    !
    ! Start time of calculation is returned to "time". 
    ! By default, "StartTimeUnit"in "NAMELIST#timeset_nml" is units of 
    ! the time. 
    ! If a unit is specified to an optional argument "unit", 
    ! a time according to the unit is returned. 
    !

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

    ! 日付および時刻の取り扱い
    ! Date and time handler
    !
    use dc_date_types, only: DC_DIFFTIME
                              ! 日時の差を表現するデータ型. 
                              ! Data type for difference about date and time
    use dc_date, only: DCDiffTimeCreate, EvalByUnit

    ! 宣言文 ; Declaration statements
    !
    implicit none
    real(DP), intent(out):: time
                              ! 計算開始時刻. 
                              ! Start time of calculation
    character(*), intent(in), optional:: unit
                              ! 計算開始時刻の単位. 
                              ! Unit of start time of calculation

    ! 作業変数
    ! Work variables
    !
    type(DC_DIFFTIME):: dcdiff_time

    ! 実行文 ; Executable statement
    !

    if ( .not. present(unit) ) then
      time = StartTimeValue
    else
      call DCDiffTimeCreate( dcdiff_time, StartTimeValue, StartTimeUnit )    ! (in)
      time = EvalByUnit( dcdiff_time, unit )
    end if

  end subroutine TimesetGetStartTime
Subroutine :

timeset モジュールの初期化を行います. NAMELIST#timeset_nml の読み込みはこの手続きで行われます.

"timeset" module is initialized. NAMELIST#timeset_nml is loaded in this procedure.

This procedure input/output NAMELIST#timeset_nml .

[Source]

  subroutine TimesetInit
    !
    ! timeset モジュールの初期化を行います. 
    ! NAMELIST#timeset_nml の読み込みはこの手続きで行われます. 
    !
    ! "timeset" module is initialized. 
    ! NAMELIST#timeset_nml is loaded in this procedure. 
    !

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

    ! NAMELIST ファイル入力に関するユーティリティ
    ! Utilities for NAMELIST file input
    !
    use namelist_util, only: namelist_filename, NmlutilMsg

    ! ファイル入出力補助
    ! File I/O support
    !
    use dc_iounit, only: FileOpen

    ! 種別型パラメタ
    ! Kind type parameter
    !
    use dc_types, only: STDOUT ! 標準出力の装置番号. Unit number of standard output

    ! 日付および時刻の取り扱い
    ! Date and time handler
    !
    use dc_date_types, only: DC_DIFFTIME
                              ! 日時の差を表現するデータ型. 
                              ! Data type for difference about date and time
    use dc_date, only: DCDiffTimeCreate, EvalSec, EvalByUnit, mod, toChar, operator(*), operator(==), operator(<), operator(>), operator(/), operator(+), operator(-)

    ! CPU 時間計測
    ! CPU time monitor
    !
    use dc_clock, only: DCClockCreate, DCClockStart

    ! 宣言文 ; Declaration statements
    !
    implicit none

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

    type(DC_DIFFTIME):: dcdiff_start
                              ! 計算開始時刻. 
                              ! Start time of calculation
    type(DC_DIFFTIME):: dcdiff_end
                              ! 計算終了時刻. 
                              ! End time of calculation
    type(DC_DIFFTIME):: dcdiff_delt
                              ! $ \Delta t $
    type(DC_DIFFTIME):: dcdiff_predict
                              ! 終了予測日時表示間隔. 
                              ! Interval of predicted end time output

    ! NAMELIST 変数群
    ! NAMELIST group name
    !
    namelist /timeset_nml/ StartTimeValue,  StartTimeUnit, EndTimeValue,    EndTimeUnit, DelTimeValue,    DelTimeUnit, PredictIntValue, PredictIntUnit, CpuTimeMoniter
          !
          ! デフォルト値については初期化手続 "timeset#TimesetInit" 
          ! のソースコードを参照のこと. 
          !
          ! Refer to source codes in the initialization procedure
          ! "timeset#TimesetInit" for the default values. 
          !

    ! 実行文 ; Executable statement
    !

    if ( timeset_inited ) return
    call InitCheck

    ! デフォルト値の設定
    ! Default values settings
    !
    Nstep = 0
    Cstep = 0

    StartTimeValue = 0.0_DP
    StartTimeUnit = 'day'

    EndTimeValue = 7.0_DP
    EndTimeUnit = 'day'

    DelTime = 180.0_DP
    DelTimeValue = 30.0_DP
    DelTimeUnit = 'min'

    PredictIntStep = 1
    PredictIntValue = 1.0_DP
    PredictIntUnit = 'day'

    CpuTimeMoniter = .true.

    ! NAMELIST の読み込み
    ! NAMELIST is input
    !
    if ( trim(namelist_filename) /= '' ) then
      call FileOpen( unit_nml, namelist_filename, mode = 'r' ) ! (in)

      rewind( unit_nml )
      read( unit_nml, nml = timeset_nml, iostat = iostat_nml ) ! (out)
      close( unit_nml )

      call NmlutilMsg( iostat_nml, module_name ) ! (in)
      if ( iostat_nml == 0 ) write( STDOUT, nml = timeset_nml )
    end if

    ! DC_DIFFTIME 型変数の設定
    ! Configure DC_DIFFTIME type variables
    !
    call DCDiffTimeCreate( dcdiff_start, StartTimeValue, StartTimeUnit )       ! (in)
    call DCDiffTimeCreate( dcdiff_end, EndTimeValue, EndTimeUnit )           ! (in)
    call DCDiffTimeCreate( dcdiff_delt, DelTimeValue, DelTimeUnit )           ! (in)
    call DCDiffTimeCreate( dcdiff_predict, PredictIntValue, PredictIntUnit)      ! (in)

    ! 時刻の正当性のチェック
    ! Check validation of time
    !
    call TimeValidCheck( dcdiff_start, dcdiff_end, dcdiff_delt, dcdiff_predict )

    ! 総ステップ数算出
    ! Calculate number of full steps
    !
    Nstep = int( ( dcdiff_end - dcdiff_start ) / dcdiff_delt )
    Nstep = Nstep + 1

    ! $ \Delta t $ [s] 算出
    ! Calculate $ \Delta t $ [s] 
    !
    DelTime = EvalSec( dcdiff_delt )

    ! 終了予測日時表示間隔ステップ数算出
    ! Calculate interval steps of predicted end time output
    PredictIntStep = int( dcdiff_predict / dcdiff_delt )
    PredictIntStep = max( PredictIntStep, 1 )

    ! DC_DIFFTIME 型変数の設定
    ! Configure DC_DIFFTIME type variables
    !
    call DCDiffTimeCreate( DCdiffCurrentTime, StartTimeValue, StartTimeUnit )            ! (in)
    call DCDiffTimeCreate( DCdiffDelTime, DelTimeValue, DelTimeUnit )            ! (in)

    ! CPU 時間計測開始 (モデル全体)
    ! Start CPU time monitoring (for entire model)
    !
    if ( CpuTimeMoniter ) then
      call DCClockCreate( clocks(clk_proc_num + 1), 'total' ) ! (in)
      call DCClockStart( clocks(clk_proc_num + 1) ) ! (in)
      clocks_name(clk_proc_num + 1) = 'total'
      clk_proc_num = clk_proc_num + 1
    end if

    ! 印字 ; Print
    !
    call MessageNotify( 'M', module_name, '----- Initialization Messages -----' )
    call MessageNotify( 'M', module_name, '  Nstep = %d', i = (/ Nstep /) )
    call MessageNotify( 'M', module_name, '  StartTime  = %f [%c]', d = (/ StartTimeValue /), c1 = trim(StartTimeUnit) )
    call MessageNotify( 'M', module_name, '  EndTime    = %f [%c]', d = (/ EndTimeValue /), c1 = trim(EndTimeUnit) )
    call MessageNotify( 'M', module_name, '  DelTime    = %f [%c]', d = (/ DelTimeValue /), c1 = trim(DelTimeUnit) )
    call MessageNotify( 'M', module_name, '             = %f [%c]', d = (/ DelTime /), c1 = 'sec' )
    call MessageNotify( 'M', module_name, '  PredictInt = %f [%c]', d = (/ PredictIntValue /), c1 = trim(PredictIntUnit) )
    call MessageNotify( 'M', module_name, '  PredictIntStep = %d', i = (/ PredictIntStep /) )
    call MessageNotify( 'M', module_name, '  CpuTimeMoniter = %b', l = (/ CpuTimeMoniter /) )
    call MessageNotify( 'M', module_name, '-- version = %c', c1 = trim(version) )

    Cstep = 1
    timeset_inited = .true.
  end subroutine TimesetInit
Subroutine :
DelTime_Value :real(DP), intent(in)
: $ Delta t $ . 単位は DelTimeUnit にて指定. Unit is specified by "DelTimeUnit".
DelTime_Unit :character(*), intent(in)
: $ Delta t $ の単位. Unit of $ Delta t $
IntStep :integer, intent(out)
: ステップ間隔. Interval of steps

与えられた DelTime_Value と DelTime_Unit を timeset モジュール内で保存している $ Delta t $ で割り, その結果を ステップ間隔として IntStep に返します.

1 以下になる場合には 1 を返します.

Given "DelTime_Value" and "DelTime_Unit" are divided by $ Delta t $ stored in "timeset" module, and the result are returned to "IntStep" as interval of step.

If the result is less than 1, 1 is returned.

[Source]

  subroutine TimesetIntStepEval( DelTime_Value, DelTime_Unit, IntStep )
    !
    ! 与えられた DelTime_Value と DelTime_Unit を
    ! timeset モジュール内で保存している $ \Delta t $ で割り, その結果を
    ! ステップ間隔として IntStep に返します. 
    !
    ! 1 以下になる場合には 1 を返します. 
    !
    ! Given "DelTime_Value" and "DelTime_Unit" are
    ! divided by $ \Delta t $ stored in "timeset" module, 
    ! and the result are returned to "IntStep" as interval of step.
    !
    ! If the result is less than 1, 1 is returned.
    !

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

    ! 日付および時刻の取り扱い
    ! Date and time handler
    !
    use dc_date_types, only: DC_DIFFTIME
                              ! 日時の差を表現するデータ型. 
                              ! Data type for difference about date and time
    use dc_date, only: DCDiffTimeCreate, EvalSec, EvalByUnit, mod, toChar, operator(*), operator(==), operator(<), operator(>), operator(/), operator(+), operator(-)

    ! 宣言文 ; Declaration statements
    !
    implicit none
    real(DP), intent(in):: DelTime_Value
                              ! $ \Delta t $ .  単位は DelTimeUnit にて指定.
                              ! Unit is specified by "DelTimeUnit". 
    character(*), intent(in):: DelTime_Unit
                              ! $ \Delta t $ の単位. 
                              ! Unit of $ \Delta t $ 
    integer, intent(out):: IntStep
                              ! ステップ間隔. 
                              ! Interval of steps

    ! 作業変数
    ! Work variables
    !
    type(DC_DIFFTIME):: dcdiff_delt_ext
                              ! $ \Delta t $
    type(DC_DIFFTIME):: dcdiff_delt_int
                              ! $ \Delta t $

    ! 実行文 ; Executable statement
    !

    ! DC_DIFFTIME 型変数の設定
    ! Configure DC_DIFFTIME type variables
    !
    call DCDiffTimeCreate( dcdiff_delt_ext, DelTime_Value, DelTime_Unit )          ! (in)
    call DCDiffTimeCreate( dcdiff_delt_int, DelTimeValue, DelTimeUnit )            ! (in)

    ! IntStep の算出
    ! Calculate IntStep
    !
    IntStep = int( dcdiff_delt_ext / dcdiff_delt_int )
    IntStep = max( IntStep, 1 )

  end subroutine TimesetIntStepEval
Subroutine :

timeset モジュール内の時刻を進めます. また, TimesetProgress#PredictIntStep で設定された値に応じて, 現在までの計算時間と計算終了予測時刻を表示します.

Progress time configured in "timeset" module. And, computation time until now and predicted end of computation time are printed according to configured "TimesetProgress#PredictIntStep"

[Source]

  subroutine TimesetProgress
    !
    ! timeset モジュール内の時刻を進めます. 
    ! また, TimesetProgress#PredictIntStep で設定された値に応じて, 
    ! 現在までの計算時間と計算終了予測時刻を表示します. 
    !
    ! Progress time configured in "timeset" module. 
    ! And, computation time until now and 
    ! predicted end of computation time are printed 
    ! according to configured "TimesetProgress#PredictIntStep"
    !

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

    ! CPU 時間計測
    ! CPU time monitor
    !
    use dc_clock, only: DCClockPredict, DCClockStop, DCClockClose, operator(+), operator(-)

    ! 日付および時刻の取り扱い
    ! Date and time handler
    !
    use dc_date, only: operator(==), operator(+)

    ! 宣言文 ; Declaration statements
    !
    implicit none

    ! 作業変数
    ! Work variables
    !
    type(CLOCK):: clock_tmp

    ! 実行文 ; Executable statement
    !

    ! 終了予測日時表示
    ! Print predicted end time
    !
    if ( mod( Cstep, PredictIntStep ) == 0 ) then
      if ( CpuTimeMoniter ) then
        clock_tmp = clocks(1)
        call DCClockStop( clock_tmp ) ! (in)
        call DCClockPredict( clock_tmp, real(Cstep) / real(Nstep) ) ! (in)
        call DCClockClose( clock_tmp ) ! (in)
      else
        call MessageNotify( 'M', module_name, ' Current/Total = [ %d / %d ]', i = (/ Cstep, Nstep /) )
      end if
    end if

    ! 時刻の進行
    ! Progress time
    !
    Cstep = Cstep + 1
    DCdiffCurrentTime = DCdiffCurrentTime + DCdiffDelTime
  end subroutine TimesetProgress
timeset_inited
Variable :
timeset_inited = .false. :logical, save, public
: 初期設定フラグ. Initialization flag

Private Instance methods

CpuTimeMoniter
Variable :
CpuTimeMoniter :logical, save
: CPU 時間計測のオンオフ On/off of CPU time monitoring
DCdiffCurrentTime
Variable :
DCdiffCurrentTime :type(DC_DIFFTIME), save
: 現在の日時. Current date and time
DCdiffDelTime
Variable :
DCdiffDelTime :type(DC_DIFFTIME), save
: $ Delta t $
DelTimeUnit
Variable :
DelTimeUnit :character(TOKEN), save
: $ Delta t $ の単位. Unit of $ Delta t $
DelTimeValue
Variable :
DelTimeValue :real(DP), save
: $ Delta t $ . 単位は DelTimeUnit にて指定. Unit is specified by "DelTimeUnit".
EndTimeUnit
Variable :
EndTimeUnit :character(TOKEN), save
: 計算開始時刻の単位. Unit of end time of calculation
EndTimeValue
Variable :
EndTimeValue :real(DP), save
: 計算終了時刻. End time of calculation
Subroutine :

依存モジュールの初期化チェック

Check initialization of dependency modules

[Source]

  subroutine InitCheck
    !
    ! 依存モジュールの初期化チェック
    !
    ! Check initialization of dependency modules

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

    ! NAMELIST ファイル入力に関するユーティリティ
    ! Utilities for NAMELIST file input
    !
    use namelist_util, only: namelist_util_inited

    ! 実行文 ; Executable statement
    !

    if ( .not. namelist_util_inited ) call MessageNotify( 'E', module_name, '"namelist_util" module is not initialized.' )

  end subroutine InitCheck
PredictIntStep
Variable :
PredictIntStep :integer, save
: 終了予測日時表示間隔ステップ数. Interval steps of predicted end time output
PredictIntUnit
Variable :
PredictIntUnit :character(TOKEN), save
: 終了予測日時表示間隔 (単位). Unit for interval of predicted end time output
PredictIntValue
Variable :
PredictIntValue :real(DP), save
: 終了予測日時表示間隔. Interval of predicted end time output
StartTimeUnit
Variable :
StartTimeUnit :character(TOKEN), save
: 計算開始時刻の単位. Unit of start time of calculation
StartTimeValue
Variable :
StartTimeValue :real(DP), save
: 計算開始時刻. Start time of calculation
Subroutine :
dcdiff_start :type(DC_DIFFTIME)
: 計算開始時刻. Start time of calculation
dcdiff_end :type(DC_DIFFTIME)
: 計算終了時刻. End time of calculation
dcdiff_delt :type(DC_DIFFTIME)
: $ Delta t $ [s]
dcdiff_predict :type(DC_DIFFTIME)
: 終了予測日時表示間隔. Interval of predicted end time output

時刻情報についての有効性をチェックします.

Check validation about infomation time

[Source]

  subroutine TimeValidCheck( dcdiff_start, dcdiff_end, dcdiff_delt, dcdiff_predict )
    !
    ! 時刻情報についての有効性をチェックします. 
    !
    ! Check validation about infomation time
    !

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

    ! 日付および時刻の取り扱い
    ! Date and time handler
    !
    use dc_date_types, only: DC_DIFFTIME
                              ! 日時の差を表現するデータ型. 
                              ! Data type for difference about date and time
    use dc_date, only: DCDiffTimeCreate, EvalSec, EvalByUnit, mod, toChar, operator(*), operator(==), operator(<), operator(>), operator(/), operator(+), operator(-)

    implicit none

    ! 作業変数
    ! Work variables
    !
    type(DC_DIFFTIME):: dcdiff_start
                              ! 計算開始時刻. 
                              ! Start time of calculation
    type(DC_DIFFTIME):: dcdiff_end
                              ! 計算終了時刻. 
                              ! End time of calculation
    type(DC_DIFFTIME):: dcdiff_delt
                              ! $ \Delta t $ [s]
    type(DC_DIFFTIME):: dcdiff_predict
                              ! 終了予測日時表示間隔. 
                              ! Interval of predicted end time output

    ! 実行文 ; Executable statement
    !

    if ( .not. 0.0_DP < EvalSec( dcdiff_end - dcdiff_start )  ) then
      call MessageNotify( 'E', module_name, 'StartTime=<%f[%c]> is later than EndTime=<%f[%c]>', d = (/ StartTimeValue, EndTimeValue /), c1 = trim(StartTimeUnit), c2 = trim(EndTimeUnit) )
    end if

    if ( dcdiff_delt > ( dcdiff_end - dcdiff_start ) ) then
      call MessageNotify( 'E', module_name, 'DelTime=<%f[%c]> is larger than ' // 'EndTime=<%f[%c]> - StartTime=<%f[%c]>.', d = (/ DelTimeValue, EndTimeValue, StartTimeValue /), c1 = trim(DelTimeUnit), c2 = trim(EndTimeUnit), c3 = trim(StartTimeUnit) )
    end if

    if ( .not. EvalSec( dcdiff_delt ) > 0.0_DP ) then
      call MessageNotify( 'E', module_name, 'DelTime=<%f[%c]> must be more than 0.', d = (/ DelTimeValue /), c1 = trim(DelTimeUnit) )
    end if

    if ( .not. EvalSec( dcdiff_predict ) > 0.0_DP ) then
      call MessageNotify( 'E', module_name, 'PredictInt=<%f[%c]> must be more than 0.', d = (/ PredictIntValue /), c1 = trim(PredictIntUnit) )
    end if

  end subroutine TimeValidCheck
clk_proc_num
Variable :
clk_proc_num = 0 :integer, save
: CPU 時間計測を行っているプロセスの数. Number of processes monitored CPU time
clkmax
Constant :
clkmax = 64 :integer, parameter
: CPU 時間計測を行うプロセスの最大数. Maximum number of processes monitored CPU time
clocks
Variable :
clocks(1:clkmax) :type(CLOCK), save
: CPU 時間計測用構造体 Derived type for monitoring CPU time
clocks_name
Variable :
clocks_name(1:clkmax) = ’’ :character(TOKEN), save
: CPU 時間計測を行っているプロセスの名称 Names of processes monitored CPU time
module_name
Constant :
module_name = ‘timeset :character(*), parameter
: モジュールの名称. Module name
version
Constant :
version = ’$Name: dcpam5-20080731 $’ // ’$Id: timeset.f90,v 1.1.1.1 2008-07-30 08:41:33 morikawa Exp $’ :character(*), parameter
: モジュールのバージョン Module version

[Validate]