Class | initialdata_baker1998 |
In: |
../src/env/initialdata_baker1998.f90
|
Baker et al. (1998) を模した初期場を作成するためのモジュール
Subroutine : | |||
z_Temp(kmin:kmax) : | real(DP), intent(out)
| ||
z_Press(kmin:kmax) : | real(DP), intent(out)
|
* Baker and Shubert (1989) で使われた初期値を再現する. * 大気安定度から温度・圧力の基本場を作る * 温位のデータは安定度で示されているので, それを温位で焼き直す * 温位から, 基本場の温度・圧力を計算する.
subroutine initialdata_baker1998_basic( z_Temp, z_Press ) ! !== 概要 ! * Baker and Shubert (1989) で使われた初期値を再現する. ! * 大気安定度から温度・圧力の基本場を作る ! * 温位のデータは安定度で示されているので, それを温位で焼き直す ! * 温位から, 基本場の温度・圧力を計算する. ! implicit none real(DP), intent(out) :: z_Press(kmin:kmax)!圧力 real(DP), intent(out) :: z_Temp(kmin:kmax) !温度 real(DP) :: r_DPTempDz(kmin:kmax) real(DP) :: z_PTemp(kmin:kmax) real(DP), parameter :: r_dPTempDz_45km = 3.9d-3 real(DP), parameter :: r_dPTempDz_48km = 0.0d0 real(DP), parameter :: r_dPTempDz_55km = 0.0d0 real(DP), parameter :: r_dPTempDz_60km = 8.35d-3 integer :: k ! 初期化 ! z_Press = 0.0d0 z_Temp = 0.0d0 ! 大気安定度を与える. ! BS1998 の表 1 を目で読み取り, それを直線近似している. ! do k = kmin, kmax if ( r_Z(k) < 45.0d3 ) then r_DPTempDz(k) = r_DPTempDz_45km + 0.30d-6 * ( r_Z(k) - 45.0d3 ) elseif ( 45.0d3 <= r_Z(k) .AND. r_Z(k) < 48.0d3 ) then r_DPTempDz(k) = r_DPTempDz_45km - 1.30d-6 * ( r_Z(k) - 45.0d3 ) elseif ( 48.0d3 <= r_Z(k) .AND. r_Z(k) < 55.0d3 ) then r_DPTempDz(k) = r_DPTempDz_48km elseif ( 55.0d3 <= r_Z(k) .AND. r_Z(k) < 60.0d3 ) then r_DPTempDz(k) = r_DPTempDz_55km + 1.67d-6 * ( r_Z(k) - 55.0d3 ) elseif ( 60.0d3 <= r_Z(k) ) then r_DPTempDz(k) = r_DPTempDz_60km - 1.20d-6 * ( r_Z(k) - 60.0d3 ) end if end do ! 温位勾配 DPTempDz を積分し温度場 z_PTemp を計算 ! 傾きは半整数格子点の値を用いる. ! z_PTemp(nz) = TempTop - r_DPTempDz(1) * dz * 0.5d0 do k = nz-1, 1, -1 z_PTemp(k) = z_PTemp(k+1) - r_DPTempDz(k) * dz end do !---------------------------------------------- ! 圧力と温度の計算 ! z_Press(nz) = 2.20d4 z_Temp(nz) = 268.0d0 ! z_Press(nz) = PressTop + (Grav * PressTop * dz * 5.0d-1) / (GasRDry * TempTop) ! z_Temp(nz) = z_PTemp(nz) * (z_Press(nz) / PressBasis) ** (GasRDry / CpDry) do k = nz-1, 1, -1 z_Press(k) = z_Press(k+1) + (Grav * z_Press(k+1) * dz) / (GasRDry * z_Temp(k+1)) z_Temp(k) = z_PTemp(k) * ((z_Press(k) / PressBasis) ** (GasRDry / CpDry)) end do end subroutine Initialdata_baker1998_basic