Class CFLCheck
In: ../src/utils/cflcheck.f90

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

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

Methods

Included Modules

dc_types dc_message gridset axesset timeset

Public Instance methods

Subroutine :
pyz_VelX(imin:imax,jmin:jmax,kmin:kmax) :real(DP), intent(in)

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

[Source]

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

    !暗黙の型宣言禁止
    implicit none
    
    !内部変数
    real(DP), intent(in) :: pyz_VelX(imin:imax,jmin:jmax,kmin:kmax)
    real(DP)             :: CFL

    !CFL 条件を求める
    CFL = (2.0d0 * DelTimeLong) * maxval(abs(pyz_VelX/xyz_dx))
  
    !メッセージ出力
    call MessageNotify( "M", module_name, "Courant number of VelX for DelTimeLong = %f", d=(/CFL/) )

!    if (CFL > 1.0d0) stop 

  end subroutine CFLCheckTimeLongVelX
Subroutine :
xqz_VelY(imin:imax,jmin:jmax,kmin:kmax) :real(DP), intent(in)

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

[Source]

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

    !暗黙の型宣言禁止
    implicit none
    
    !内部変数
    real(DP), intent(in) :: xqz_VelY(imin:imax,jmin:jmax,kmin:kmax)
    real(DP)             :: CFL

    !CFL 条件を求める
    CFL = (2.0d0 * DelTimeLong) * maxval(abs(xqz_VelY/xyz_dy))
  
    !メッセージ出力
    call MessageNotify( "M", module_name, "Courant number of VelY for DelTimeLong = %f", d=(/CFL/) )

!    if (CFL > 1.0d0) stop 

  end subroutine CFLCheckTimeLongVelY
Subroutine :
xyr_VelZ(imin:imax,jmin:jmax,kmin:kmax) :real(DP), intent(in)

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

[Source]

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

    !暗黙の型宣言禁止
    implicit none
    
    !内部変数
    real(DP), intent(in) :: xyr_VelZ(imin:imax,jmin:jmax,kmin:kmax)
    real(DP)             :: CFL
    
    !CFL 条件を求める
    CFL = (2.0d0 * DelTimeLong) * maxval(abs(xyr_VelZ/xyz_dz)) 
    
    !メッセージ出力
    call MessageNotify( "M", module_name, "Courant number of VelZ for DelTimeLong = %f", d=(/CFL/) )
    
!    if (CFL > 1.0d0) stop 

  end subroutine CFLCheckTimeLongVelZ
Subroutine :
xyz_VelSound(imin:imax,jmin:jmax, kmin:kmax) :real(DP), intent(in)
: 音速

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

[Source]

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

    !暗黙の型宣言禁止
    implicit none
    
    !変数定義
    real(DP), intent(in) :: xyz_VelSound(imin:imax,jmin:jmax, kmin:kmax)
                                        !音速
    real(DP)             :: xyz_CFL(imin:imax,jmin:jmax, kmin:kmax)
                                        !クーラン数
!    real(DP)             :: DxyzMin      !最小格子間隔
    
    !音速と CFL 条件を求める
!    DxyzMin = min(min(minval(xyz_dx), minval(xyz_dy)), minval(xyz_dz))
!    CFL = DelTimeShort * maxval(xyz_VelSound) / DxyzMin

    xyz_CFL = DelTimeShort * xyz_VelSound * ((1.0d0 / (xyz_dx ** 2.0d0) + 1.0d0 / (xyz_dy ** 2.0d0)) ** 0.5d0)

    !メッセージ
    call MessageNotify( "M", module_name, "Sound Wave Velocity = %f", d=(/maxval(xyz_VelSound)/) )
!    call MessageNotify( "M", &
!      & module_name, &
!      & "min(DelX, DelY, DelZ) = %f", d=(/DxyzMin/) )
    call MessageNotify( "M", module_name, "DelTimeShort = %f", d=(/DelTimeShort/) )

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