Class timeset
In: ../src/setup/timeset.f90

引数に与えられた NAMELIST ファイルから, 時刻に関する情報を取得し, 保管するための変数型モジュール

Methods

Included Modules

dc_types dc_iounit dc_message mpi_wrapper namelist_util

Public Instance methods

DelTimeLong
Variable :
DelTimeLong = 2.0d0 :real(DP), save, public
: 長いタイムステップ
DelTimeOutput
Variable :
DelTimeOutput = 2.0d0 :real(DP), save, public
: 出力タイムステップ
DelTimeShort
Variable :
DelTimeShort = 2.0d-1 :real(DP), save, public
: 短いタイムステップ
EndTime
Variable :
EndTime = 3600.0d0 :real(DP), save, public
: 計算終了時刻
FlagInitialRun
Variable :
FlagInitialRun = .false. :logical, save, public
IntegPeriod
Variable :
IntegPeriod = 3600.0d0 :real(DP), save, public
: 積分時間
NstepOutput
Variable :
NstepOutput = 20 :integer, save, public
: リスタートファイルへの出力
NstepShort
Variable :
NstepShort = 20 :integer, save, public
: 短いタイムステップのステップ数
RestartTime
Variable :
RestartTime = 0.0d0 :real(DP), save, public
: 計算開始時刻
TimeA
Variable :
TimeA :real(DP), save, public
: 時刻 t + del t
TimeB
Variable :
TimeB :real(DP), save, public
: 時刻 t - del t
TimeN
Variable :
TimeN :real(DP), save, public
: 時刻 t
Subroutine :

[Source]

  subroutine TimesetDelTimeHalf    

    implicit none

    ! 最初の一回目はオイラースキームで回すので, 
    ! 長い時間刻みを半分にし, 短い時間で回すループ回数も半分にする.
    ! Asselin の時間フィルタの係数をゼロとする. 
    !
    TimeB       = TimeN
    DelTimeLong = DelTimeLongSave * 0.5d0
    NstepShort  = NstepShortSave  / 2
    tfil        = 0.0d0

    !---------------------------------------------------------------
    ! 確認
    !
    if (myrank == 0) then
      call MessageNotify( "M", "timeset_init", "Initial DelTimeLong  = %f", d=(/DelTimeLong/) )
      call MessageNotify( "M", "timeset_init", "Initial NstepShort   = %d", i=(/NstepShort/) )
!      call MessageNotify( "M", &
!        & "timeset_init", "Asselin Time Filter coefficient = %d", d=(/tfil/) )
!      write(*,*) tfil
    end if    

  end subroutine TimesetDelTimeHalf
Subroutine :

[Source]

  subroutine TimesetProgress

    implicit none

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

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

    ! Asselin の時間フィルタの係数
    tfil = tfilSave

  end subroutine TimesetProgress
tfil
Variable :
tfil = 1.0d-1 :real(DP), save, public
: アセリンの時間フィルタの係数
Subroutine :

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

This procedure input/output NAMELIST#timeset_nml .

[Source]

  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
    
    ! 時刻・タイムステップの設定
    !   実数の割り算なので, 念の為に. 
    !    
    NstepShort = 2 * ( nint( DelTimeLong * 1.0d2 ) / nint( DelTimeShort * 1.0d2 ) )

    ! 時間刻みの設定
    !   元々の値を保管しておく
    DelTimeLongSave = DelTimeLong
    NstepShortSave  = NstepShort

    ! リスタートファイルを書き出すタイミング
    !
    NstepOutput = nint( DelTimeOutput * 1.0d2 ) / nint( DelTimeLong * 1.0d2 )

    ! 時刻の設定
    !
    TimeB = RestartTime - DelTimeLong
    TimeN = RestartTime
    TimeA = RestartTime + DelTimeLong      
    
    ! リスタートか否か. リスタートなら  .false.
    !
    if ( nint( RestartTime * 1.0d2 ) == 0 ) then 
      FlagInitialRun = .true.
    else
      FlagInitialRun = .false.
    end if

    !---------------------------------------------------------------
    ! 確認
    !
    if (myrank == 0) then
      call MessageNotify( "M", "timeset_init", "DelTimeLong  = %f", d=(/DelTimeLong/) )
      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