| Class | CFLCheck |
| In: |
util/cflcheck.f90
|
CFL 条件のチェックをするためのパッケージ型モジュール
* 音波に対して CFL 条件をチェック * 入力された速度に対して CFL 条件をチェック
| Subroutine : | |
| pz_VelX(DimXMin:DimXMax, DimZMin:DimZMax) : | real(8), intent(in) |
水平速度に対して CFL 条件をチェック.
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 : | |
| xr_VelZ(DimXMin:DimXMax, DimZMin:DimZMax) : | real(8), intent(in) |
水平速度に対して CFL 条件をチェック.
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
| Subroutine : | |||
| xz_VelSound(DimXMin:DimXMax, DimZMin:DimZMax) : | real(8), intent(in)
|
音波に対して CFL 条件をチェック
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) / DelX
!メッセージ
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