!---------------------------------------------------------------------
!     Copyright (C) GFD Dennou Club, 2004, 2006. All rights reserved.
!---------------------------------------------------------------------
!= Subroutine CFLCheck
!
!   * Developer: SUGIYAMA Ko-ichiro
!   * Version: $Id: cflcheck.f90,v 1.2 2005/04/22 15:01:52 sugiyama Exp $ 
!   * Tag Name: $Name:  $
!   * Change History: 
!
!== Overview 
!
!CFL 条件のチェックをするためのパッケージ型モジュール
!  * 音波に対して CFL 条件をチェック
!  * 入力された速度に対して CFL 条件をチェック
!
!== Error Handling
!
!== Known Bugs
!
!== Note
!
!== Future Plans
!
!Error Handling を gt4f90io を利用するように変更
!

module CFLCheck
  !
  !CFL 条件のチェックをするためのパッケージ型モジュール
  !  * 音波に対して CFL 条件をチェック
  !  * 入力された速度に対して CFL 条件をチェック
  !

  !モジュール読み込み
  use gridset, only: DimXMin,       &! x 方向の配列の下限
    &                DimXMax,       &! x 方向の配列の上限
    &                DimZMin,       &! z 方向の配列の下限
    &                DimZMax,       &! z 方向の配列の上限
    &                DelX,          &! X 方向の格子点間隔
    &                DelZ            ! Z 方向の格子点間隔
  use timeset, only: DelTimeShort,  &!短い時間ステップ
    &                DelTimeLong     !長い時間ステップ
  
  !暗黙の型宣言禁止
  implicit none

  !private 属性の指定
  private

  !関数を public 属性に設定
  public CFLCheckTimeShort
  public CFLCheckTimeLongVelX
  public CFLCheckTimeLongVelZ
  
contains  

!!!-----------------------------------------------------------------!!!
  subroutine CFLCheckTimeShort( xz_VelSound )
    !
    !音波に対して CFL 条件をチェック
    ! 

    !暗黙の型宣言禁止
    implicit none
    
    !変数定義
    real(8), intent(in) :: xz_VelSound(DimXMin:DimXMax, DimZMin:DimZMax)
                                        !音速
    real(8)             :: CFL          !クーラン数
    
    !音速と CFL 条件を求める
    CFL = DelTimeShort / min(DelX / maxval(xz_VelSound), DelZ / maxval(xz_VelSound))

    !メッセージ
    write(*,*) "CFL Condition for DelTimeSort: ", CFL

    !警告メッセージ
    if ( CFL >= 1.0) then 
      write(*,*) "CFL Condition is broken, DelTimeShort * VelSound > Del"
      write(*,*) "sound wave velocity ", maxval(xz_VelSound)
      write(*,*) "Del                 ", min(DelX, DelZ)
      write(*,*) "DelTimeShort        ", DelTimeShort
      stop
    end if
  
  end subroutine CFLCheckTimeShort


!!!-----------------------------------------------------------------!!!
  subroutine CFLCheckTimeLongVelX( pz_VelX )
    !
    !水平速度に対して CFL 条件をチェック. 
    ! 

    !暗黙の型宣言禁止
    implicit none
    
    !内部変数
    real(8), intent(in) :: pz_VelX(DimXMin:DimXMax, DimZMin:DimZMax)
    real(8)             :: CFL
    
    !音速と CFL 条件を求める
    CFL = DelTimeLong/ max(DelX / maxval(pz_VelX), - DelX / minval(pz_VelX)) 
  
    !メッセージ出力
    write(*,*) "CFL Condition of VelX for DelTimeLong: ", CFL

  end subroutine CFLCheckTimeLongVelX
    

!!!-----------------------------------------------------------------!!!
  subroutine CFLCheckTimeLongVelZ( xr_VelZ )
    !
    !水平速度に対して CFL 条件をチェック. 
    ! 

    !暗黙の型宣言禁止
    implicit none
    
    !内部変数
    real(8), intent(in) :: xr_VelZ(DimXMin:DimXMax, DimZMin:DimZMax)
    real(8)             :: CFL
    
    !音速と CFL 条件を求める
    CFL = DelTimeLong / max(DelZ / maxval(xr_VelZ), - DelZ / minval(xr_VelZ)) 
    
    !メッセージ出力
    write(*,*) "CFL Condition of VelZ for DelTimeLong: ", CFL
    
  end subroutine CFLCheckTimeLongVelZ
  
end module CFLCheck
