!= 時刻管理
!
!= Time control
!
! Authors::   Yasuhiro MORIKAWA
! Version::   $Id: gridset.f90,v 1.5 2007/06/09 13:12:30 sugiyama Exp $ 
! Tag Name::  $Name:  $
! Copyright:: Copyright (C) GFD Dennou Club, 2008. All rights reserved.
! License::   See COPYRIGHT[link:../../../COPYRIGHT]
!

module timeset
  !
  != 時刻管理
  !
  != Time control
  !
  ! <b>Note that Japanese and English are described in parallel.</b>
  !
  ! 時刻情報の管理を行う. 
  !
  ! Information of time is controlled. 
  !

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

  ! 種別型パラメタ
  ! Kind type parameter
  !
  use dc_types, only: DP, &  ! 倍精度実数型. Double precision. 
    &                 TOKEN  ! キーワード.   Keywords. 

  ! メッセージ出力
  ! Message output
  !
  use dc_message, only: MessageNotify

  ! 宣言文 ; Declaration statements
  !
  implicit none
  private

  ! 公開手続き
  ! Public procedure
  !
  public:: TimesetInit, TimesetClose, TimesetProgress
  public:: TimesetIntStepEval
  public:: TimesetClockStart, TimesetClockStop

  ! 公開変数
  ! Public variables
  !
  logical, save, public:: initialized = .false.
                              ! 初期設定フラグ. 
                              ! Initialization flag
  integer, save, public:: Nstep = 0
                              ! 総ステップ数. Number of full steps
  integer, save, public:: Cstep = 0
                              ! 現在のステップ数. Current steps

  real(DP), save, public:: StartTimeValue = 0.0_DP
                              ! 計算開始時刻. 
                              ! Start time of calculation
  character(TOKEN), save, public:: StartTimeUnit = 'min'
                              ! 計算開始時刻の単位. 
                              ! Unit of start time of calculation

  real(DP), save, public:: EndTimeValue = 7.0_DP
                              ! 計算終了時刻. 
                              ! End time of calculation
  character(TOKEN), save, public:: EndTimeUnit = 'day'
                              ! 計算開始時刻の単位. 
                              ! Unit of end time of calculation

  real(DP), save, public:: DelTime = 180.0_DP
                              ! $ \Delta t $ [s]
  real(DP), save, public:: DelTimeValue = 30.0_DP
                              ! $ \Delta t $ .  単位は DelTimeUnit にて指定.
                              ! Unit is specified by "DelTimeUnit". 
  character(TOKEN), save, public:: DelTimeUnit = 'min'
                              ! $ \Delta t $ の単位. 
                              ! Unit of $ \Delta t $ 

  integer, save, public:: PredictIntStep = 1
                              ! 終了予測日時表示間隔ステップ数. 
                              ! Interval steps of predicted end time output
  real(DP), save, public:: PredictIntValue = 1.0_DP
                              ! 終了予測日時表示間隔. 
                              ! Interval of predicted end time output
  character(TOKEN), save, public:: PredictIntUnit = 'day'
                              ! 終了予測日時表示間隔 (単位). 
                              ! Unit for interval of predicted end time output

  ! 非公開変数
  ! Private variables
  !

  namelist /timeset_nml/ &
    & StartTimeValue,  StartTimeUnit, &
    & EndTimeValue,    EndTimeUnit, &
    & DelTimeValue,    DelTimeUnit, &
    & PredictIntValue, PredictIntUnit

  character(*), parameter:: module_name = 'timeset'
                              ! モジュールの名称. 
                              ! Module name
  character(*), parameter:: version = &
    & '$Name:  $' // &
    & '$Id:  $'
                              ! モジュールのバージョン
                              ! Module version

  ! INTERFACE 文 ; INTERFACE statements
  !
  interface TimesetInit
    module procedure TimesetInit
  end interface

  interface TimesetClose
    module procedure TimesetClose
  end interface

  interface TimesetProgress
    module procedure TimesetProgress
  end interface

  interface TimesetIntStepEval
    module procedure TimesetIntStepEval
  end interface

  interface TimesetClockStart
    module procedure TimesetClockStart
  end interface

  interface TimesetClockStop
    module procedure TimesetClockStop
  end interface

contains

  subroutine TimesetInit
    !
    ! timeset モジュールの初期化を行います. 
    ! NAMELIST#timeset_nml の値が読み込まれます. 
    !
    ! "timeset" module is initialized. 
    ! Values of NAMELIST#timeset_nml are loaded. 
    !

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

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

    ! 実行文 ; Executable statement
    !

    if ( initialized ) return
    call InitCheck

    ! NAMELIST の読み込み
    ! NAMELIST is input
    !
    call FileOpen( unit_nml, &          ! (out)
      & namelist_filename, mode = 'r' ) ! (in)

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

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

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

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

    ! 総ステップ数算出
    ! 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_delt / dcdiff_predict )
    PredictIntStep = max( PredictIntStep, 1 )

    ! 印字 ; 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 /) )

    Cstep = 1
    initialized = .true.
  end subroutine TimesetInit

  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
    !

    ! 宣言文 ; Declaration statements
    !
    implicit none

    ! 作業変数
    ! Work variables
    !

    ! 実行文 ; Executable statement
    !

    ! 終了予測日時表示
    ! Print predicted end time
    !
    if ( mod( Cstep, PredictIntStep ) == 0 ) then
      write(*,*) 'it is planned that prediction time is printed'
    end if

    ! 時刻の進行
    ! Progress time
    !
    Cstep = Cstep + 1

  end subroutine TimesetProgress

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

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

    ! 宣言文 ; Declaration statements
    !
    implicit none

    ! 作業変数
    ! Work variables
    !

    ! 実行文 ; Executable statement
    !

    write(*,*) 'it is planned that total computation time is printed'

  end subroutine TimesetClose

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

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

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

    ! 作業変数
    ! Work variables
    !

    ! 実行文 ; Executable statement
    !

  end subroutine TimesetClockStart

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

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

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

    ! 作業変数
    ! Work variables
    !

    ! 実行文 ; Executable statement
    !

  end subroutine TimesetClockStop

  subroutine TimesetIntStepEval( &
    & DelTime_Value, DelTime_Unit, & ! (in)
    & IntStep &                      ! (out)
    & )
    !
    ! 与えられた 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
    !

    ! 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(-)

    ! 宣言文 ; 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, &  ! (out)
      & DelTime_Value, DelTime_Unit )          ! (in)
    call DCDiffTimeCreate( dcdiff_delt_int, &  ! (out)
      & DelTimeValue, DelTimeUnit )            ! (in)

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

  end subroutine TimesetIntStepEval

  subroutine TimeValidCheck( &
    & dcdiff_start, dcdiff_end, dcdiff_delt, dcdiff_predict & ! (in)
    & )
    !
    ! 時刻情報についての有効性をチェックします. 
    !
    ! 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, StartTimeValue, EndTimeValue /), &
        & c1 = trim(DelTimeUnit), &
        & c2 = trim(StartTimeUnit), c3 = trim(EndTimeUnit) )
    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

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

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

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

    ! 実行文 ; Executable statement
    !

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

  end subroutine InitCheck

end module timeset
