Class CFLCheck
In: util/cflcheck.f90

CFL 条件のチェックをするためのパッケージ型モジュール

 * 音波に対して CFL 条件をチェック
 * 入力された速度に対して CFL 条件をチェック

Methods

Included Modules

dc_message debugset gridset timeset

Public Instance methods

Subroutine :
pz_VelX(DimXMin:DimXMax, DimZMin:DimZMax) :real(8), intent(in)

水平速度に対して CFL 条件をチェック.

[Source]

  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)) 
  
    !メッセージ出力
    call MessageNotify( "M", "CFLCheckTimeLongVelX", "Courant number of VelX for DelTimeLong = %f", d=(/CFL/) )

  end subroutine CFLCheckTimeLongVelX
Subroutine :
xr_VelZ(DimXMin:DimXMax, DimZMin:DimZMax) :real(8), intent(in)

水平速度に対して CFL 条件をチェック.

[Source]

  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)) 
    
    !メッセージ出力
    call MessageNotify( "M", "CFLCheckTimeLongVelZ", "Courant number of VelZ for DelTimeLong = %f", d=(/CFL/) )
  end subroutine CFLCheckTimeLongVelZ
Subroutine :
xz_VelSound(DimXMin:DimXMax, DimZMin:DimZMax) :real(8), intent(in)
: 音速

音波に対して CFL 条件をチェック

[Source]

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

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

    !メッセージ
    if (cpurank == 0) then 
      call MessageNotify( "M", "CFLCheckTimeShort", "Sound Wave Velocity = %f", d=(/maxval(xz_VelSound)/) )
      call MessageNotify( "M", "CFLCheckTimeShort", "min(DelX, DelZ) = %f", d=(/min(DelX, DelZ)/) )
      call MessageNotify( "M", "CFLCheckTimeShort", "DelTimeShort = %f", d=(/DelTimeSHort/) )
    end if

    !警告メッセージ
    if (cpurank == 0) then 
      if ( CFL >= 1.0) then 
        call MessageNotify( "E", "CFLCheckTimeShort", "CFL Condition is broken, DelTimeShort * VelSound > min(DelX, DelZ)")
      else
        call MessageNotify( "M", "CFLCheckTimeShort", "Courant number for DelTimeSort = %f", d=(/CFL/) )
      end if
    end if
  
  end subroutine CFLCheckTimeShort