| Class | initialdata_sounding |
| In: |
../src/env/initialdata_sounding.f90
|
サウンディングファイルから初期場を決めるためのサブルーチン サウンディングファイルはテキストファイルで書かれている. 将来的には netCDF に変更する予定.
| Subroutine : | |||
| z_Temp(kmin:kmax) : | real(DP), intent(out)
| ||
| z_Press(kmin:kmax) : | real(DP), intent(out)
|
subroutine initialdata_sounding_basic( z_Temp, z_Press )
!
!
implicit none
real(DP), intent(out):: z_Press(kmin:kmax) !圧力
real(DP), intent(out):: z_Temp(kmin:kmax) !温度
real(DP) :: r_Press(kmin:kmax) !圧力
real(DP) :: r_Temp(kmin:kmax) !温度
real(DP) :: z_DTempDZ(kmin:kmax)
real(DP) :: DTempDZ
integer :: i, k, k1, k2
logical :: flag
! 初期化
!
z_Temp = 0.0d0
r_Temp = 0.0d0
z_Press = 0.0d0
r_Press = 0.0d0
flag = .false.
! データ読み込み
!
do k = kmin, kmax
do i = 1, NumRec
if ( r_Z(k) == r_tmpAlt(i) ) then
r_Temp(k) = r_tmpTemp(i)
r_Press(k) = r_tmpPress(i)
end if
end do
end do
! データ読み込み
! うまく半格子に合う場合もあるので, その場合はデータを優先.
!
do k = kmin, kmax
do i = 1, NumRec
if ( z_Z(k) == r_tmpAlt(i) ) then
flag = .true.
z_Temp(k) = r_tmpTemp(i)
z_Press(k) = r_tmpPress(i)
end if
end do
end do
! r => z の変換
!
if (.NOT. flag) then
do k = kmin+1, kmax
z_Temp(k) = ( r_Temp(k-1) + r_Temp(k) ) / 2.0d0
z_Press(k) = ( r_Press(k-1) + r_Press(k) ) / 2.0d0
end do
end if
! 確認
!
! do k = kmin, kmax
! write(*,*) z_z(k), z_Temp(k), z_Press(k)
! end do
!!!
!!! 圏界面が存在する場合には, そこでの温度勾配を緩やかにする.
!!!
! 現在の温度分布の勾配を取る
!
do k = 1, kmax
z_DTempDZ(k) = (z_Temp(k) - z_Temp(k-1)) / dz
end do
! 対流圏界面より上の扱い. 指定された高度の温度減率を使い続ける.
!
k1 = minloc( z_Z, 1, z_Z > AltTr )
! 温度減率
!
DTempDZ = z_DTempDZ(k1-1)
k2 = DelAlt / dz
do k = k1, kmax
! 温度源率をゼロに近づける.
!
DTempDZ = min( -1.0d-14, DTempDZ - DTempDZ / k2 * ( k - k1 ) )
!基本場の温度を決める
z_Temp(k) = z_Temp(k-1) + DTempDZ * dz
!圧力を静水圧平衡から計算
z_Press(k) = z_Press(k-1) * ( ( z_Temp(k-1) / z_Temp(k) ) ** (Grav / ( DTempDZ * GasRDry ) ) )
end do
end subroutine Initialdata_sounding_basic
| Subroutine : |
設定ファイルから出力ファイルに記載する情報を読み込む
This procedure input/output NAMELIST#initialdata_sounding_nml .
subroutine initialdata_sounding_init
!
!設定ファイルから出力ファイルに記載する情報を読み込む
!
!暗黙の型宣言禁止
implicit none
!内部変数
integer :: AltCol = 0
integer :: TempCol = 0
integer :: PressCol = 0
integer :: VelXCol = 0
integer :: VelYCol = 0
integer :: unit !設定ファイル用装置番号
integer :: io
character(30) :: SoundingFile
integer, parameter :: maxch=12
character(len=100) :: buf, eachcol(maxch)
integer :: num, MaxCol
!設定ファイルから読み込む出力ファイル情報
!
NAMELIST /initialdata_sounding_nml/ SoundingFile, AltCol, TempCol, PressCol, VelXCol, VelYCol, AltTr, DelAlt
!設定ファイルから出力ファイルに記載する情報を読み込む
!
call FileOpen(unit, file=namelist_filename, mode='r')
read(unit, NML=initialdata_sounding_nml)
close(unit)
! write(*,*) SoundingFile, AltCol, TempCol, PressCol
! 初期化
!
io = 0
NumRec = 0
MaxCol = max( AltCol, max( TempCol, PressCol ) )
! write(*,*) "MaxCol", MaxCol
r_tmpAlt = 0.0d0
r_tmpTemp = 0.0d0
r_tmpPress = 0.0d0
r_tmpVelX = 0.0d0
r_tmpVelY = 0.0d0
! ファイルのオープン
!
open (17, file=SoundingFile, status='old')
! ファイル呼び出し
!
do while ( io == 0 )
! 1 行分読み出し
!
read (17, '(a)', IOSTAT=io) buf
! 行をカンマ区切りで分割
!
call devidecsv( buf, eachcol, maxch, num )
! 確認
!
! write(*,*) num
! do i=1, num
! write(*,*) i, eachcol(i)(1:len_trim(eachcol(i)))
! end do
! num の値が小さいものはヘッダとみなす.
!
if (num >= MaxCol) then
! 行数の計算
!
NumRec = NumRec + 1
! 値の代入
!
if (AltCol > 0) read( eachcol(AltCol)(1:len_trim(eachcol(AltCol))), *) r_tmpAlt(NumRec)
if (TempCol > 0) read( eachcol(TempCol)(1:len_trim(eachcol(TempCol))), *) r_tmpTemp(NumRec)
if (PressCol > 0) read( eachcol(PressCol)(1:len_trim(eachcol(PressCol))), *) r_tmpPress(NumRec)
if (VelXCol > 0) read( eachcol(VelXCol)(1:len_trim(eachcol(VelXCol))), *) r_tmpVelX(NumRec)
if (VelYCol > 0) read( eachcol(VelYCol)(1:len_trim(eachcol(VelYCol))), *) r_tmpVelY(NumRec)
end if
end do
! 確認
!
! write(*,*) z_tmpAlt(1:NumRec)
! ファイルのクローズ
!
close (17)
end subroutine initialdata_sounding_init
| Subroutine : | |
| pyz_VelX(imin:imax,jmin:jmax,kmin:kmax) : | real(DP), intent(out) |
| xqz_VelY(imin:imax,jmin:jmax,kmin:kmax) : | real(DP), intent(out) |
subroutine initialdata_sounding_wind(pyz_VelX, xqz_VelY)
implicit none
real(DP), intent(out) :: pyz_VelX(imin:imax,jmin:jmax,kmin:kmax)
real(DP), intent(out) :: xqz_VelY(imin:imax,jmin:jmax,kmin:kmax)
real(DP) :: pyr_VelX(imin:imax,jmin:jmax,kmin:kmax)
real(DP) :: xqr_VelY(imin:imax,jmin:jmax,kmin:kmax)
integer :: i, k
logical :: flag
!初期化
pyz_VelX = 0.0d0
pyr_VelX = 0.0d0
xqz_VelY = 0.0d0
xqr_VelY = 0.0d0
flag = .false.
! データ読み込み
!
do k = kmin, kmax
do i = 1, NumRec
if ( r_Z(k) == r_tmpAlt(i) ) then
pyr_VelX(:,:,k) = r_tmpVelX(i)
xqr_VelY(:,:,k) = r_tmpVelY(i)
end if
end do
end do
! データ読み込み
! うまく半格子に合う場合もあるので, その場合はデータを優先.
!
do k = kmin, kmax
do i = 1, NumRec
if ( z_Z(k) == r_tmpAlt(i) ) then
flag = .true.
pyz_VelX(:,:,k) = r_tmpVelX(i)
xqz_VelY(:,:,k) = r_tmpVelY(i)
end if
end do
end do
! r => z の変換
!
if (.NOT. flag) then
do k = kmin+1, kmax
pyz_VelX(:,:,k) = ( pyr_VelX(:,:,k-1) + pyr_VelX(:,:,k) ) * 5.0d-1
xqz_VelY(:,:,k) = ( xqr_VelY(:,:,k-1) + xqr_VelY(:,:,k) ) * 5.0d-1
end do
end if
! 確認
!
! do k = kmin, kmax
! write(*,*) z_z(k), pyz_VelX(1,1,k), xqz_VelY(1,1,k)
! end do
end subroutine initialdata_sounding_wind