Class timeset
In: setup/timeset.f90

時刻管理

Time control

Note that Japanese and English are described in parallel.

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

Information of time is controlled.

Variables List

TimeB :ステップ $ t - Delta t $ の時刻
TimeN :ステップ $ t $ の時刻
TimeA :ステップ $ t + Delta t $ の時刻
DelTime :$ Delta t $ [s] (実数型の数値)
StartTime :計算開始時刻
EndTime :計算終了時刻
StartDate :計算開始日時
StartDateValid :計算開始日時の有効性
———— :————
TimeB :Time of step $ t - Delta t $
TimeN :Time of step $ t $
TimeA :Time of step $ t + Delta t $
DelTime :$ Delta t $ [s] (numerical value of float type)
StartTime :Start time of calculation
EndTime :End time of calculation
StartDate :Start date of calculation
StartDateValid :Validation of start date of calculation

Procedures List

TimesetInit :timeset モジュールの初期化
TimesetDelTimeHalf :Δt を一時的に半分に
TimesetProgress :時刻の進行
TimesetClockStart :計算時間計測開始
TimesetClockStop :計算時間計一時停止
TimesetClose :timeset モジュールの終了処理
———— :————
TimesetInit :Initialize "timeset" module
TimesetDelTimeHalf :Reduce delta t to half temporarily
TimesetProgress :Progress time
TimesetClockStart :Start measurement of computation time
TimesetClockStop :Pause measurement of computation time
TimesetClose :Terminate "timeset" module

NAMELIST

NAMELIST#timeset_nml

Methods

Included Modules

dc_types dc_message dc_date_types dc_clock namelist_util dc_iounit dc_date dc_string

Public Instance methods

DelTime
Variable :
DelTime :real(DP), save, public
: $ Delta t $ [s]
EndTime
Variable :
EndTime :type(DC_DIFFTIME), save, public
: 計算終了時刻. End time of calculation
StartDate
Variable :
StartDate :type(DC_DATETIME), save, public
: 計算開始日時. Start date of calculation
StartDateValid
Variable :
StartDateValid :logical, save, public
: 計算開始日時の有効性. NAMELIST#timeset_nmlDate に指定が 行われた場合には .true. になります.

Validation of start date of calculation If "Date" of "NAMELIST#timeset_nml" is specified explicitly, this value becomes ".true.".

StartTime
Variable :
StartTime :type(DC_DIFFTIME), save, public
: 計算開始時刻. Start time of calculation
TimeA
Variable :
TimeA :type(DC_DIFFTIME), save, public
: ステップ $ t + Delta t $ の時刻. Time of step $ t + Delta t $.
TimeB
Variable :
TimeB :type(DC_DIFFTIME), save, public
: ステップ $ t - Delta t $ の時刻. Time of step $ t - Delta t $.
TimeN
Variable :
TimeN :type(DC_DIFFTIME), save, public
: ステップ $ t $ の時刻. Time of step $ t $.
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 :

計算の初回だけはオイラー法を用いるため, 一時的に Δt を半分にします. TimesetProgress が呼ばれた段階で Δt は元に戻ります.

Delta t is reduced to half temporarily in order to use Euler method at initial step. Delta t is returned to default, when "TimesetProgress" is called.

[Source]

  subroutine TimesetDelTimeHalf
    ! 
    ! 計算の初回だけはオイラー法を用いるため, 
    ! 一時的に Δt を半分にします. 
    ! TimesetProgress が呼ばれた段階で Δt は元に戻ります. 
    ! 
    ! Delta t is reduced to half temporarily 
    ! in order to use Euler method at initial step. 
    ! Delta t is returned to default, when "TimesetProgress" is called. 
    !
    
    ! 宣言文 ; Declaration statements
    !
    implicit none

    ! 作業変数
    ! Work variables
    !

    ! 実行文 ; Executable statement
    !
    if ( flag_half ) return
    DelTime = DelTime / 2.0_DP
    flag_half = .true.

  end subroutine TimesetDelTimeHalf
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: DCDateTimeCreate, DCDiffTimeCreate, EvalSec, EvalByUnit, mod, toChar, toCharCal, operator(*), operator(==), operator(<), operator(>), operator(/), operator(+), operator(-)

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

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

    ! 宣言文 ; 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

    real, parameter:: date_zero(8) = 0.0
                              ! NAMELIST による Date の有効性を調べるための変数
                              ! Variable for check of validation of "Date" by NAMELIST

    ! NAMELIST 変数群
    ! NAMELIST group name
    !
    namelist /timeset_nml/ StartTimeValue,  StartTimeUnit, EndTimeValue,    EndTimeUnit, DelTimeValue,    DelTimeUnit, PredictIntValue, PredictIntUnit, CpuTimeMoniter, Date, Calendar
          !
          ! デフォルト値については初期化手続 "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
    !
    StartTimeValue = 0.0_DP
    StartTimeUnit = 'day'

    EndTimeValue = 7.0_DP
    EndTimeUnit = 'day'

    DelTime = 180.0_DP
    DelTimeValue = 30.0_DP
    DelTimeUnit = 'min'
    DelTimeSave = DelTime
    flag_half = .false.

    PredictIntValue = 1.0_DP
    PredictIntUnit = 'day'
    call DCDiffTimeCreate( PredictIntTime, PredictIntValue, PredictIntUnit )     ! (in)
    call DCDiffTimeCreate( PredictPrevTime, StartTimeValue, StartTimeUnit )       ! (in)

    CpuTimeMoniter = .true.

    Date = 0.0
    Calendar = ''
    StartDateValid = .false.

    ! 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( StartTime, StartTimeValue, StartTimeUnit )       ! (in)
    call DCDiffTimeCreate( EndTime, EndTimeValue, EndTimeUnit )           ! (in)
    call DCDiffTimeCreate( DCdiffDelTime, DelTimeValue, DelTimeUnit )            ! (in)
    call DCDiffTimeCreate( PredictIntTime, PredictIntValue, PredictIntUnit )    ! (in)
    PredictPrevTime = StartTime - DCdiffDelTime

    ! 時刻の正当性のチェック
    ! Check validation of time
    !
    call TimeValidCheck( StartTime, EndTime, DCdiffDelTime, PredictIntTime )

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

    ! 各タイムステップの時刻の設定
    ! Configure time on each time step
    !
    call DCDiffTimeCreate( TimeN, StartTimeValue, StartTimeUnit )   ! (in)
    TimeB = TimeN - DCdiffDelTime
    TimeA = TimeN + DCdiffDelTime


    ! 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

    ! 計算開始日時
    ! Start date of calculation
    !
    if ( any( date_zero(1:8) /= Date(1:8) ) ) then
      call DCDateTimeCreate( StartDate, year = Date(1), mon = Date(2), day = Date(3), hour = Date(4), min = Date(5), sec = real( Date(6), DP ) , zone_hour = Date(7), zone_min = Date(8), caltype_str = Calendar )                                     ! (in)
      StartDateValid = .true.
    else
      StartDateValid = .false.
    end if

    ! 印字 ; Print
    !
    call MessageNotify( 'M', module_name, '----- Initialization Messages -----' )
    call MessageNotify( 'M', module_name, '  StartTime  = %f [%c]', d = (/ StartTimeValue /), c1 = trim(StartTimeUnit) )
    if ( StartDateValid ) then
      call MessageNotify( 'M', module_name, '  StartDate  = %c', c1 = trim(toChar(StartDate)) )
      call MessageNotify( 'M', module_name, '  Calendar   = %c', c1 = trim(toCharCal(StartDate)) )
    end if
    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, '  CpuTimeMoniter = %b', l = (/ CpuTimeMoniter /) )
    call MessageNotify( 'M', module_name, '-- version = %c', c1 = trim(version) )

    timeset_inited = .true.
  end subroutine TimesetInit
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(+), operator(/), operator(-), operator(>=), toChar

    ! 宣言文 ; Declaration statements
    !
    implicit none

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

    ! 実行文 ; Executable statement
    !

    ! Δt を元に戻す. 
    ! Delta t is returned to default
    !
    if ( flag_half ) then
      DelTime = DelTimeSave
      flag_half = .false.
    end if

    ! 終了予測日時表示
    ! Print predicted end time
    !
    if ( TimeN - PredictPrevTime >= PredictIntTime ) then
      PredictPrevTime = TimeN
      if ( CpuTimeMoniter ) then
        clock_tmp = clocks(1)
        call DCClockStop( clock_tmp ) ! (in)
        call DCClockPredict( clock_tmp, real( ( TimeA - StartTime ) / ( EndTime - StartTime ) ) ) ! (in)
        call DCClockClose( clock_tmp ) ! (in)
      else
        call MessageNotify( 'M', module_name, ' Current/Total = [ %c / %c ]', c1 = trim(toChar(TimeA - StartTime)), c2 = trim(toChar(EndTime - StartTime)) )
      end if
    end if

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

Private Instance methods

Calendar
Variable :
Calendar :character(TOKEN), save
: 暦法. Calendar
CpuTimeMoniter
Variable :
CpuTimeMoniter :logical, save
: CPU 時間計測のオンオフ On/off of CPU time monitoring
DCdiffDelTime
Variable :
DCdiffDelTime :type(DC_DIFFTIME), save
: $ Delta t $
Date
Variable :
Date(8) :integer, save
: 計算開始日時. (年月日時分秒, タイムゾーン時分) Start date of calculation. (year, month, day, hour, minute, second, and hour, minute of time zone)
DelTimeSave
Variable :
DelTimeSave :real(DP), save
: $ Delta t $ [s] のデフォルト値. ("TimesetDelTimeHalf" で使用される)

Default value of $ Delta t $ [s]. (for "TimesetDelTimeHalf")

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
PredictIntTime
Variable :
PredictIntTime :type(DC_DIFFTIME), save
: 終了予測日時表示時間間隔. Interval time 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
PredictPrevTime
Variable :
PredictPrevTime :type(DC_DIFFTIME), save
: 前回の終了予測日時表示時間. Time when predicted end time is output previously
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), intent(in)
: 計算開始時刻. Start time of calculation
dcdiff_end :type(DC_DIFFTIME), intent(in)
: 計算終了時刻. End time of calculation
dcdiff_delt :type(DC_DIFFTIME), intent(in)
: $ Delta t $ [s]
dcdiff_predict :type(DC_DIFFTIME), intent(in)
: 終了予測日時表示間隔. 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
    type(DC_DIFFTIME), intent(in):: dcdiff_start
                              ! 計算開始時刻. 
                              ! Start time of calculation
    type(DC_DIFFTIME), intent(in):: dcdiff_end
                              ! 計算終了時刻. 
                              ! End time of calculation
    type(DC_DIFFTIME), intent(in):: dcdiff_delt
                              ! $ \Delta t $ [s]
    type(DC_DIFFTIME), intent(in):: dcdiff_predict
                              ! 終了予測日時表示間隔. 
                              ! Interval of predicted end time output

    ! 作業変数
    ! Work variables
    !

    ! 実行文 ; 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
flag_half
Variable :
flag_half :logical, save
: TimesetDelTimeHalf によって $ Delta t $ が 半分になっていることを示すフラグ.

Flag that shows $ Delta t $ is reduced to half by "TimesetDelTimeHalf"

module_name
Constant :
module_name = ‘timeset :character(*), parameter
: モジュールの名称. Module name
version
Constant :
version = ’$Name: dcpam5-20081118 $’ // ’$Id: timeset.f90,v 1.7 2008-10-07 09:20:05 morikawa Exp $’ :character(*), parameter
: モジュールのバージョン Module version

[Validate]