!= Module CFLCheck_3D ! ! Authors:: SUGIYAMA Ko-ichiro ! Version:: $Id: cflcheck_3d.f90,v 1.3 2007-10-01 09:01:13 odakker Exp $ ! Tag Name:: $Name: arare4-20080627 $ ! Copyright:: Copyright (C) GFD Dennou Club, 2006. All rights reserved. ! License:: See COPYRIGHT[link:../../COPYRIGHT] ! !== Overview ! !CFL 条件のチェックをするためのパッケージ型モジュール ! * 音波に対して CFL 条件をチェック ! * 入力された速度に対して CFL 条件をチェック ! !== Error Handling ! !== Known Bugs ! !== Note ! !== Future Plans ! !Error Handling を gt4f90io を利用するように変更 ! module CFLCheck_3d ! !CFL 条件のチェックをするためのパッケージ型モジュール ! * 音波に対して CFL 条件をチェック ! * 入力された速度に対して CFL 条件をチェック ! !モジュール読み込み use dc_types, only: DP use dc_message, only: MessageNotify use gridset_3d, only: DimXMin, &! x 方向の配列の下限 & DimXMax, &! x 方向の配列の上限 & DimYMin, &! z 方向の配列の下限 & DimYMax, &! z 方向の配列の上限 & DimZMin, &! z 方向の配列の下限 & DimZMax, &! z 方向の配列の上限 & xyz_dX, &! x 方向の格子点間隔 & xyz_dY, &! y 方向の格子点間隔 & xyz_dZ ! z 方向の格子点間隔 use timeset, only: DelTimeShort, &!短い時間ステップ & DelTimeLong !長い時間ステップ !暗黙の型宣言禁止 implicit none !private 属性の指定 private !関数を public 属性に設定 public CFLCheckTimeShort public CFLCheckTimeLongVelX public CFLCheckTimeLongVelY public CFLCheckTimeLongVelZ contains !!!-----------------------------------------------------------------!!! subroutine CFLCheckTimeShort( xyz_VelSound ) ! !音波に対して CFL 条件をチェック ! !暗黙の型宣言禁止 implicit none !変数定義 real(DP), intent(in) :: xyz_VelSound & & (DimXMin:DimXMax,DimYMin:DImYMax, DimZMin:DimZMax) !音速 real(DP) :: CFL !クーラン数 real(DP) :: DxyzMin !最小格子間隔 !音速と CFL 条件を求める DxyzMin = min(min(minval(xyz_dx), minval(xyz_dy)), minval(xyz_dz)) CFL = DelTimeShort * maxval(xyz_VelSound) / DxyzMin !メッセージ call MessageNotify( "M", & & "CFLCheckTimeShort", & & "Sound Wave Velocity = %f", d=(/maxval(xyz_VelSound)/) ) call MessageNotify( "M", & & "CFLCheckTimeShort", & & "min(DelX, DelY, DelZ) = %f", d=(/DxyzMin/) ) call MessageNotify( "M", & & "CFLCheckTimeShort", & & "DelTimeShort = %f", d=(/DelTimeSHort/) ) !警告メッセージ 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 subroutine CFLCheckTimeShort !!!-----------------------------------------------------------------!!! subroutine CFLCheckTimeLongVelX( pyz_VelX ) ! !水平速度に対して CFL 条件をチェック. ! !暗黙の型宣言禁止 implicit none !内部変数 real(DP), intent(in) :: pyz_VelX & & (DimXMin:DimXMax,DimYMin:DimYMax,DimZMin:DimZMax) real(DP) :: CFL !CFL 条件を求める CFL = DelTimeLong * maxval(abs(pyz_VelX/xyz_dx)) !メッセージ出力 call MessageNotify( "M", & & "CFLCheckTimeLongVelX", & & "Courant number of VelX for DelTimeLong = %f", d=(/CFL/) ) end subroutine CFLCheckTimeLongVelX !!!-----------------------------------------------------------------!!! subroutine CFLCheckTimeLongVelY( xqz_VelY ) ! !水平速度に対して CFL 条件をチェック. ! !暗黙の型宣言禁止 implicit none !内部変数 real(DP), intent(in) :: xqz_VelY & & (DimXMin:DimXMax,DimYMin:DimYMax,DimZMin:DimZMax) real(DP) :: CFL !CFL 条件を求める CFL = DelTimeLong * maxval(abs(xqz_VelY/xyz_dy)) !メッセージ出力 call MessageNotify( "M", & & "CFLCheckTimeLongVelY", & & "Courant number of VelY for DelTimeLong = %f", d=(/CFL/) ) end subroutine CFLCheckTimeLongVelY !!!-----------------------------------------------------------------!!! subroutine CFLCheckTimeLongVelZ( xyr_VelZ ) ! !水平速度に対して CFL 条件をチェック. ! !暗黙の型宣言禁止 implicit none !内部変数 real(DP), intent(in) :: xyr_VelZ & & (DimXMin:DimXMax,DimYMin:DImYMax,DimZMin:DimZMax) real(DP) :: CFL !CFL 条件を求める CFL = DelTimeLong * maxval(abs(xyr_VelZ/xyz_dz)) !メッセージ出力 call MessageNotify( "M", & & "CFLCheckTimeLongVelZ", & & "Courant number of VelZ for DelTimeLong = %f", d=(/CFL/) ) end subroutine CFLCheckTimeLongVelZ end module CFLCheck_3d