!= Module TimeSet
!
! Authors::   SUGIYAMA Ko-ichiro, ODAKA Masatsugu
! Version::   $Id: timeset.f90,v 1.7 2011-07-15 08:43:29 sugiyama Exp $ 
! Tag Name::  $Name: arare5-20111010 $
! Copyright:: Copyright (C) GFD Dennou Club, 2006. All rights reserved.
! License::   See COPYRIGHT[link:../../COPYRIGHT]
!
!== Overview 
!
!引数に与えられた NAMELIST ファイルから, 時刻に関する情報を取得し, 
!保管するための変数型モジュール
!
!== Error Handling
!
!== Known Bugs
!
!== Note
!
!== Future Plans
!
!

module timeset
  !
  !引数に与えられた NAMELIST ファイルから, 時刻に関する情報を取得し, 
  !保管するための変数型モジュール
  !
  
  !モジュール読み込み
  use dc_types,   only: DP, STRING
  use dc_iounit,  only: FileOpen
  use dc_message, only: MessageNotify
  use mpi_wrapper, only: myrank
  use namelist_util, only: namelist_filename
  
  !暗黙の型宣言禁止
  implicit none

  !属性
  private
  
  ! Public Interface
  real(DP), save, public  :: TimeA                   !時刻 t + \del t
  real(DP), save, public  :: TimeN                   !時刻 t
  real(DP), save, public  :: TimeB                   !時刻 t - \del t
  real(DP), save, public  :: DelTimeLong  = 2.0d0    !長いタイムステップ
  real(DP), save, public  :: DelTimeLongSave  = 2.0d0 !長いタイムステップ
  real(DP), save, public  :: DelTimeShort = 2.0d-1   !短いタイムステップ
  real(DP), save, public  :: RestartTime  = 0.0d0    !計算開始時刻
  real(DP), save, public  :: IntegPeriod  = 3600.0d0 !積分時間
  real(DP), save, public  :: EndTime      = 3600.0d0 !計算終了時刻
  real(DP), save, public  :: DelTimeOutput= 2.0d0    !出力タイムステップ
  integer,  save, public  :: NstepShort = 20         !短いタイムステップのステップ数
  integer,  save, public  :: NstepShortSave = 20     !短いタイムステップのステップ数
  integer,  save, public  :: NstepOutput    = 20     !リスタートファイルへの出力
  
  ! 公開要素
  public timeset_init, TimesetProgress

contains
   
  subroutine timeset_init()
    !
    !NAMELIST から必要な情報を読み取り, 時間関連の変数の設定を行う. 
    !

    !暗黙の型宣言禁止
    implicit none

    !内部変数
    integer    :: unit

    !---------------------------------------------------------------    
    ! NAMELIST から情報を取得
    !
    NAMELIST /timeset_nml/ &
      & DelTimeLong, DelTimeShort, IntegPeriod, RestartTime, DelTimeOutput
    
    call FileOpen(unit, file=namelist_filename, mode='r')
    read(unit, NML=timeset_nml)
    close(unit)

    ! 計算終了時刻
    !
    EndTime = RestartTime + IntegPeriod
    
    ! 時刻・タイムステップの設定
    !   t=0 の時は, 最初の 1 ループだけオイラー法で解くので細工する. 
    !    
    NstepShort = 2 * nint( DelTimeLong / DelTimeShort )

    DelTimeLongSave = DelTimeLong
    NstepShortSave  = NstepShort

    NstepOutput = nint( DelTimeOutput / DelTimeLong )

    if (RestartTime /= 0.0d0) then 
      TimeB = RestartTime - DelTimeLong
      TimeN = RestartTime
      TimeA = RestartTime + DelTimeLong      
    else
      TimeB = RestartTime
      TimeN = RestartTime
      TimeA = RestartTime + DelTimeLong
      DelTimeLong = DelTimeLong * 5.0d-1
      NstepShort = NstepShort / 2
    end if

    !---------------------------------------------------------------
    ! 確認
    !
    if (myrank == 0) then
    
      !長い時間ステップが短い時間ステップで割り切れない場合には警告を出す
      if(mod(DelTimeLong, DelTimeShort) /= 0) then 
        call MessageNotify( "W", &
          & "timeset_init", "mod(DelTimeLong, DelTimeShort) is not zero")
      end if

      !長い時間ステップが短い時間ステップで割り切れない場合には警告を出す
      if(mod(DelTimeOutput, DelTimeLong) /= 0) then 
        call MessageNotify( "W", &
          & "timeset_init", "mod(DelTimeOutput, DelTimeLong) is not zero")
      end if
      
      call MessageNotify( "M", &
        & "timeset_init", "DelTimeLong  = %f", d=(/DelTimeLongSave/) )
      call MessageNotify( "M", &
        & "timeset_init", "DelTimeShort = %f", d=(/DelTimeShort/) )
      call MessageNotify( "M", &
        & "timeset_init", "Restarttime  = %f", d=(/Restarttime/)  )
      call MessageNotify( "M", &
        & "timeset_init", "IntegPeriod  = %f", d=(/IntegPeriod/) )
      call MessageNotify( "M", &
        & "timeset_init", "EndTime      = %f", d=(/EndTime/) )
      call MessageNotify( "M", &
        & "timeset_init", "DelTimeOutput= %f", d=(/DelTimeOutput/) )
      call MessageNotify( "M", &
        & "timeset_init", "NstepShort   = %d", i=(/NstepShort/) )
      call MessageNotify( "M", &
        & "timeset_init", "NstepOutput  = %d", i=(/NstepOutput/) )
    end if
    
  end subroutine timeset_init


  subroutine TimesetProgress

    implicit none

    ! 時刻刻み幅を直す 
    !
    DelTimeLong = DelTimeLongSave
    NstepShort  = NstepShortSave

    ! 時刻を進める
    !
    TimeB = TimeN
    TimeN = TimeA
    TimeA = TimeA + DelTimeLong

  end subroutine TimesetProgress
  
end module timeset
