| Class | Thermo_Advanced_Function |
| In: |
thermo_advanced_function.f90
|
基礎ルーチン, 関数集を use して複雑な熱力学関数を計算するモジュール
| Function : | |||
| CAPE : | real | ||
| tmp_p(:) : | real, intent(in)
| ||
| tmp_z(size(tmp_p)) : | real, intent(in)
| ||
| tmp_qv(size(tmp_p)) : | real, intent(in)
| ||
| tmp_t(size(tmp_p)) : | real, intent(in)
| ||
| z_ref : | real, intent(in)
| ||
| undeff : | real, intent(in), optional |
real function CAPE( tmp_p, tmp_z, tmp_qv, tmp_t, z_ref, undeff )
use Thermo_Const
use Thermo_Function
use Algebra
use Statistics
implicit none
real, intent(in) :: tmp_p(:) ! 気圧 [Pa]
real, intent(in) :: tmp_z(size(tmp_p)) ! 高度 [m]
real, intent(in) :: tmp_qv(size(tmp_p)) ! 混合比 [kg/kg]
real, intent(in) :: tmp_t(size(tmp_p)) ! 温度 [K]
real, intent(in) :: z_ref ! パーセルを持ち上げる基準高度 [m]
real, intent(in), optional :: undeff
integer :: nx, i, j, iz_LFC, iz_LNB
real :: PLFC, TLFC, PLNB, ZLFC, ZLNB
real :: tmp(size(tmp_p)), t_par(size(tmp_p))
t_par=0.0
nx=size(tmp_p)
!-- LFC, LNB での高度, 温度, 圧力の計算
TLFC=T_LFC( z_ref, tmp_z, tmp_t, tmp_p, tmp_qv )
ZLFC=z_LFC( z_ref, tmp_z, tmp_t, tmp_p, tmp_qv )
ZLNB=z_LNB( z_ref, tmp_z, tmp_t, tmp_p, tmp_qv )
call interpo_search_1d( tmp_z, ZLFC, iz_LFC, undeff )
call interpo_search_1d( tmp_z, ZLNB, iz_LNB, undeff )
call interpolation_1d( tmp_z(iz_LFC), tmp_z(iz_LFC+1), tmp_p(iz_LFC), tmp_p(iz_LFC+1), ZLFC, PLFC )
call interpolation_1d( tmp_z(iz_LNB), tmp_z(iz_LNB+1), tmp_p(iz_LNB), tmp_p(iz_LNB+1), ZLNB, PLNB )
do i=1,nx
if(tmp_p(i)<=PLFC.and.tmp_p(i)>=PLNB)then
if(present(undeff))then
if(undeff==tmp_t(i))then
tmp(i)=undeff
else
t_par(i)=moist_laps_temp( PLFC, TLFC, tmp_p(i) )
tmp(i)=(t_par(i)-tmp_t(i))/tmp_p(i)
end if
else
t_par(i)=moist_laps_temp( PLFC, TLFC, tmp_p(i) )
tmp(i)=(t_par(i)-tmp_t(i))/tmp_p(i)
end if
end if
end do
call rectangle_int( tmp_p, tmp, PLNB, PLFC, CAPE, undeff )
! PLNB が先に来ているが, rectangle_int の処理の都合でそうなっている.
! 積分の向きは, \int^{PLNB}_{PLFC} で行われる.
CAPE=-CAPE*Rd
return
end function
| Function : | |||
| CIN : | real | ||
| tmp_p(:) : | real, intent(in)
| ||
| tmp_z(size(tmp_p)) : | real, intent(in)
| ||
| tmp_qv(size(tmp_p)) : | real, intent(in)
| ||
| tmp_t(size(tmp_p)) : | real, intent(in)
| ||
| z_ref : | real, intent(in)
| ||
| undeff : | real, intent(in), optional |
real function CIN( tmp_p, tmp_z, tmp_qv, tmp_t, z_ref, undeff )
use Thermo_Const
use Algebra
use Thermo_Function
use Statistics
implicit none
real, intent(in) :: tmp_p(:) ! 気圧 [Pa]
real, intent(in) :: tmp_z(size(tmp_p)) ! 高度 [m]
real, intent(in) :: tmp_qv(size(tmp_p)) ! 混合比 [kg/kg]
real, intent(in) :: tmp_t(size(tmp_p)) ! 温度 [K]
real, intent(in) :: z_ref ! パーセルを持ち上げる基準高度 [m]
real, intent(in), optional :: undeff
integer :: nx, i, j, iz_LFC, iz_ref, iz_LCL
real :: PLFC, TLFC, ZLFC, PLCL, ZLCL, TLCL, p_ref, T_ref, qv_ref
real :: tmp(size(tmp_p)), t_par(size(tmp_p))
t_par=0.0
nx=size(tmp_p)
!-- z_ref を基準とした, LCL, LFC での温度, 高度, 圧力の計算
!-- LFC 高度の変数計算
TLFC=T_LFC( z_ref, tmp_z, tmp_t, tmp_p, tmp_qv )
ZLFC=z_LFC( z_ref, tmp_z, tmp_t, tmp_p, tmp_qv )
call interpo_search_1d( tmp_z, ZLFC, iz_LFC, undeff )
call interpolation_1d( tmp_z(iz_LFC), tmp_z(iz_LFC+1), tmp_p(iz_LFC), tmp_p(iz_LFC+1), ZLFC, PLFC )
!-- reference 高度の変数計算
call interpo_search_1d( tmp_z, z_ref, iz_ref, undeff )
call interpolation_1d( tmp_z(iz_ref), tmp_z(iz_ref+1), tmp_p(iz_ref), tmp_p(iz_ref+1), z_ref, p_ref )
call interpolation_1d( tmp_z(iz_ref), tmp_z(iz_ref+1), tmp_t(iz_ref), tmp_t(iz_ref+1), z_ref, T_ref )
call interpolation_1d( tmp_z(iz_ref), tmp_z(iz_ref+1), tmp_qv(iz_ref), tmp_qv(iz_ref+1), z_ref, qv_ref )
!-- LCL 高度の変数計算
TLCL=TqvP_2_TLCL( t_ref, qv_ref, p_ref )
ZLCL=z_LCL( z_ref, tmp_z, tmp_t, tmp_p, tmp_qv )
call interpo_search_1d( tmp_z, ZLCL, iz_LCL, undeff )
call interpolation_1d( tmp_z(iz_LCL), tmp_z(iz_LCL+1), tmp_p(iz_LCL), tmp_p(iz_LCL+1), ZLCL, PLCL )
! 周囲の温度は T で得られているので, まず LCL を計算し, その高度
! までを乾燥断熱で, それ以降を湿潤断熱でパーセル温度を計算.
do i=1,nx
if(tmp_z(i)>=z_ref.and.tmp_z(i)<=ZLCL)then
t_par(i)=T_ref-(g/Cpd)*(tmp_z(i)-z_ref)
else
t_par(i)=moist_laps_temp( PLCL, TLCL, tmp_p(i) )
end if
if(present(undeff))then
if(undeff==tmp_t(i))then
tmp(i)=undeff
else
tmp(i)=(t_par(i)-tmp_t(i))/tmp_p(i)
end if
else
tmp(i)=(t_par(i)-tmp_t(i))/tmp_p(i)
end if
end do
!do i=1,nx
! if(tmp_p(i)<p_ref.and.tmp_p(i)>PLFC)then
! t_par(i)=t_bot-(g/CPD)*(height(i+1)-height(i))
! if(PLCL>PLFC)then ! 強い逆転層のある場合
! if(p(i)<=PLCL.and.p(i)>=PLFC)then
! call moist_laps_calc( PLCL, TLCL, p(i), t_par(i) )
! end if
! end if
! end if
!end do
call rectangle_int( tmp_p, tmp, PLFC, p_ref, CIN, undeff )
! PLNB が先に来ているが, rectangle_int の処理の都合でそうなっている.
! 積分の向きは, \int^{PLNB}_{PLFC} で行われる.
CIN=-CIN*Rd
return
end function
| Function : | |||
| Louis : | real | ||
| z : | real, intent(in)
| ||
| z0m : | real, intent(in)
| ||
| richard : | real, intent(in)
|
Louis(1980) で提案されている大気の不安定度を考慮したバルク係数の補正係数を計算する関数
real function Louis( z, z0m, richard )
! Louis(1980) で提案されている大気の不安定度を考慮したバルク係数の補正係数を計算する関数
use Thermo_Const
implicit none
real, intent(in) :: z ! cm を求める高度 [m]
real, intent(in) :: z0m ! モデルで計算される粗度高度 [m]
real, intent(in) :: richard ! バルクリチャードソン数
real, parameter :: b=5.0, c=5.0
real :: cm_tmp, zratio
cm_tmp=(kalm/(log(z)-log(z0m)))**2
zratio=z/z0m
if(richard<0.0)then
Louis=1.0-((2.0*b*richard)/(1.0+3.0*b*c*cm_tmp*sqrt(-richard*zratio)))
else
Louis=1.0/(1.0+2.0*b*richard*sqrt(1.0+c*richard))
end if
return
end function
| Function : | |||
| Rich : | real | ||
| za : | real, intent(in)
| ||
| pta : | real, intent(in)
| ||
| ptg : | real, intent(in)
| ||
| va : | real, intent(in)
| ||
| qva : | real, intent(in)
| ||
| qvs : | real, intent(in)
|
バルクリチャードソン数を計算する関数
real function Rich( za, pta, ptg, va, qva, qvs ) ! バルクリチャードソン数を計算する関数 use Phys_Const use Thermo_Function implicit none real, intent(in) :: za ! リチャードソン数を計算する高度 [m] real, intent(in) :: pta ! za での仮温位 [K] real, intent(in) :: ptg ! 地表面での温位 [K] real, intent(in) :: va ! 高度 za での水平風速の絶対値 [m/s] real, intent(in) :: qva ! za での混合比 [kg/kg] real, intent(in) :: qvs ! 地表面での飽和混合比 [kg/kg] real :: ptvg, ptva, dpt ptvg=ptg*((1.0+eps_rdrv*qvs)/(1.0+qvs)) ptva=pta*((1.0+eps_rdrv*qva)/(1.0+qva)) dpt=ptva-ptvg Rich=(g*za*dpt)/(ptva*(va**2)) return end function
| Function : | |||
| T_LFC : | real | ||
| z_ref : | real, intent(in)
| ||
| z(:) : | real, intent(in)
| ||
| temp(size(z)) : | real, intent(in)
| ||
| pres(size(z)) : | real, intent(in)
| ||
| qv(size(z)) : | real, intent(in)
|
LFC 高度での温度を計算する. ある基準高度での相当温位を上昇させたとき, 周囲の飽和相当温位と一致する高度が LFC である. ここでは, 周囲の飽和相当温位と交差する 2 点を求め, そこから高度方向に 内挿することで高度を決定する.
real function T_LFC( z_ref, z, temp, pres, qv ) ! LFC 高度での温度を計算する.
! ある基準高度での相当温位を上昇させたとき, 周囲の飽和相当温位と一致する高度が
! LFC である.
! ここでは, 周囲の飽和相当温位と交差する 2 点を求め, そこから高度方向に
! 内挿することで高度を決定する.
use Statistics
use Thermo_Function
implicit none
real, intent(in) :: z_ref ! 基準高度 [m]
real, intent(in) :: z(:) ! 高度座標 [m]
real, intent(in) :: temp(size(z)) ! 温度 [K]
real, intent(in) :: pres(size(z)) ! 圧力 [Pa]
real, intent(in) :: qv(size(z)) ! 水蒸気混合比 [kg/kg]
integer :: i, j, k, nz, iz, iept
real :: sept(size(z))
real :: eptiz, eptiz1, ept_ref, z_sol
nz=size(z)
!-- 鉛直方向の飽和相当温位を計算.
do i=1,nz
sept(i)=thetaes_Bolton( temp(i), pres(i) )
end do
!-- z_ref での相当温位を計算する.
call interpo_search_1d( z, z_ref, iz )
!-- iz と iz+1 での相当温位を計算.
eptiz=thetae_Bolton( temp(iz), qv(iz), pres(iz) )
eptiz1=thetae_Bolton( temp(iz+1), qv(iz+1), pres(iz+1) )
!-- これらの間から相当温位を内挿
call interpolation_1d( z(iz), z(iz+1), eptiz, eptiz1, z_ref, ept_ref )
!-- z_ref より上の高度で ept_ref と飽和相当温位の関係を計算.
do i=iz+1,nz
if((sept(i)-ept_ref)*(sept(i+1)-ept_ref)<0.0)then
iept=i
exit
end if
end do
!-- パーセルの相当温位が環境場の飽和相当温位と逆転する直前の高度が iept である.
!-- ここで, iept と iept+1 の飽和相当温位から高度を内挿する.
call interpolation_1d( sept(iept), sept(iept+1), temp(iept), temp(iept+1), ept_ref, z_sol )
T_LFC=z_sol
return
end function
| Function : | |||
| T_LNB : | real | ||
| z_ref : | real, intent(in)
| ||
| z(:) : | real, intent(in)
| ||
| temp(size(z)) : | real, intent(in)
| ||
| pres(size(z)) : | real, intent(in)
| ||
| qv(size(z)) : | real, intent(in)
| ||
| opt : | integer, intent(in), optional
|
LNB 高度での温度を計算する. LFC 高度以上で再び環境場の飽和相当温位と交わる高度が LNB である. ただし, 実際の測定データでは, 測器の微小変動により, LFC の直上で LNB に 到達する場合がある. そこで, opt として, LFC を求めてから LNB を計算する方法と 上空から下向きに下ろして最初に交差する高度を LNB と定義する手法の 2 パターン用意することにする. opt = 1, LFC から求める. opt = 2 直上から計算. デフォルトでは opt = 1 を計算.
real function T_LNB( z_ref, z, temp, pres, qv, opt ) ! LNB 高度での温度を計算する.
! LFC 高度以上で再び環境場の飽和相当温位と交わる高度が LNB である.
! ただし, 実際の測定データでは, 測器の微小変動により, LFC の直上で LNB に
! 到達する場合がある.
! そこで, opt として, LFC を求めてから LNB を計算する方法と
! 上空から下向きに下ろして最初に交差する高度を LNB と定義する手法の
! 2 パターン用意することにする.
! opt = 1, LFC から求める. opt = 2 直上から計算.
! デフォルトでは opt = 1 を計算.
use Statistics
use Thermo_Function
implicit none
real, intent(in) :: z_ref ! 基準高度 [m]
real, intent(in) :: z(:) ! 高度座標 [m]
real, intent(in) :: temp(size(z)) ! 温度 [K]
real, intent(in) :: pres(size(z)) ! 圧力 [Pa]
real, intent(in) :: qv(size(z)) ! 水蒸気混合比 [kg/kg]
integer, intent(in), optional :: opt ! 計算手法のオプション
integer :: i, j, k, nz, iz, iept, counter
real :: sept(size(z))
real :: eptiz, eptiz1, ept_ref, z_sol
nz=size(z)
!-- 鉛直方向の飽和相当温位を計算.
do i=1,nz
sept(i)=thetaes_Bolton( temp(i), pres(i) )
end do
!-- z_ref での相当温位を計算する.
call interpo_search_1d( z, z_ref, iz )
!-- iz と iz+1 での相当温位を計算.
eptiz=thetae_Bolton( temp(iz), qv(iz), pres(iz) )
eptiz1=thetae_Bolton( temp(iz+1), qv(iz+1), pres(iz+1) )
!-- これらの間から相当温位を内挿
call interpolation_1d( z(iz), z(iz+1), eptiz, eptiz1, z_ref, ept_ref )
counter=0
if(present(opt))then
if(opt==2)then
do i=nz,iz+1,-1 ! 上から下に下ろす.
if((sept(i)-ept_ref)*(sept(i-1)-ept_ref)<0.0)then
iept=i-1 ! 上から下に下ろしているので, 1 要素下のデータが iept.
exit
end if
end do
else
!-- z_ref より上の高度で ept_ref と飽和相当温位の関係を計算.
do i=iz+1,nz
if((sept(i)-ept_ref)*(sept(i+1)-ept_ref)<0.0)then
counter=counter+1
if(counter==2)then
iept=i
exit
end if
end if
end do
end if
else
do i=iz+1,nz
if((sept(i)-ept_ref)*(sept(i+1)-ept_ref)<0.0)then
counter=counter+1
if(counter==2)then
iept=i
exit
end if
end if
end do
end if
!-- パーセルの相当温位が環境場の飽和相当温位と逆転する直前の高度が iept である.
!-- ここで, iept と iept+1 の飽和相当温位から高度を内挿する.
call interpolation_1d( sept(iept), sept(iept+1), temp(iept), temp(iept+1), ept_ref, z_sol )
T_LNB=z_sol
return
end function
| Function : | |||
| cm : | real | ||
| z : | real, intent(in)
| ||
| z0m : | real, intent(in)
| ||
| richard : | real, intent(in), optional
|
運動量に関するバルク係数を計算する関数
real function cm( z, z0m, richard )
! 運動量に関するバルク係数を計算する関数
use Thermo_Const
implicit none
real, intent(in) :: z ! cm を求める高度 [m]
real, intent(in) :: z0m ! モデルで計算される粗度高度 [m]
real, intent(in), optional :: richard ! Louis (1980) のスキームで計算する場合のバルクリチャードソン数
if(present(richard))then
cm=(kalm/(log(z)-log(z0m)))**2
cm=cm*Louis( z, z0m, richard )
else
cm=(kalm/(log(z)-log(z0m)))**2
end if
return
end function
| Function : | |||
| precip_water : | real | ||
| p(:) : | real, intent(in)
| ||
| qv(size(p)) : | real, intent(in)
| ||
| undef : | real, intent(in), optional
|
可降水量を計算する. 単位は [kg/m^2] 積分範囲は p の要素数の上下端で自動的に指定. 気象がわかる数と式より, 高度座標で積分するのではなく, 静力学平衡から気圧座標へ置き直して積分.
real function precip_water( p, qv, undef ) ! 可降水量を計算する. 単位は [kg/m^2]
! 積分範囲は p の要素数の上下端で自動的に指定.
! 気象がわかる数と式より, 高度座標で積分するのではなく,
! 静力学平衡から気圧座標へ置き直して積分.
use Algebra
implicit none
real, intent(in) :: p(:) ! 圧力 [Pa]
real, intent(in) :: qv(size(p)) ! 水蒸気混合比 [kg/kg]
real, intent(in), optional :: undef ! undef
integer :: nx, i
real, dimension(size(p)) :: tmp_p, tmp_qv
real :: precip
nx=size(p)
if(present(undef))then
do i=1,nx
if(p(i)==undef.or.qv(i)==undef)then
tmp_qv(i)=0.0
tmp_p(i)=tmp_p(i-1)
else
tmp_qv(i)=qv(i)
tmp_p(i)=p(i)
end if
end do
else
do i=1,nx
tmp_qv(i)=qv(i)
tmp_p(i)=p(i)
end do
end if
call rectangle_int( tmp_p, tmp_qv, tmp_p(nx), tmp_p(1), precip )
write(*,*) "precip", precip, tmp_qv(1), tmp_qv(nx)
precip_water=-precip/g
return
end function
| Function : | |||
| ust : | real | ||
| cmd : | real, intent(in)
| ||
| va : | real, intent(in)
|
摩擦速度 u_* を計算する関数
real function ust( cmd, va ) ! 摩擦速度 u_* を計算する関数 implicit none real, intent(in) :: cmd ! 高度 za でのバルク係数 real, intent(in) :: va ! 高度 za での水平風の絶対値 [m/s] ust=va*sqrt(cmd) end function
| Function : | |||
| z_LCL : | real | ||
| z_ref : | real, intent(in)
| ||
| z(:) : | real, intent(in)
| ||
| temp(size(z)) : | real, intent(in)
| ||
| pres(size(z)) : | real, intent(in)
| ||
| qv(size(z)) : | real, intent(in)
|
LCL 高度を計算する. パーセルが z_ref から断熱過程で上昇し, この温度に達したときの高度が LCL 高度となるので, z_ref から乾燥断熱減率で以下のように計算する. LCL 高度を z_LCL とし, パーセルの基準高度を z_ref とすると, LCL に達するまでは乾燥断熱減率 Gamma _d で変化するので, $Gamma _d=\frac{g}{C_p} =-frac{T_LCL-T_ref}{z_LCL-z_ref} $ という式が成り立つ. ここで, z_LCL について解くと, $z_LCL=z_ref+frac{C_p}{g} (T_ref-T_LCL)$ となる.
real function z_LCL( z_ref, z, temp, pres, qv ) ! LCL 高度を計算する.
! パーセルが z_ref から断熱過程で上昇し, この温度に達したときの高度が
! LCL 高度となるので, z_ref から乾燥断熱減率で以下のように計算する.
! LCL 高度を z_LCL とし, パーセルの基準高度を z_ref とすると,
! LCL に達するまでは乾燥断熱減率 \Gamma _d で変化するので,
! $\Gamma _d=\frac{g}{C_p} =-\frac{T_LCL-T_ref}{z_LCL-z_ref} $
! という式が成り立つ. ここで, z_LCL について解くと,
! $z_LCL=z_ref+\frac{C_p}{g} (T_ref-T_LCL)$
! となる.
use Thermo_Const
use Phys_Const
use Statistics
implicit none
real, intent(in) :: z_ref ! 基準高度 [m]
real, intent(in) :: z(:) ! 高度座標 [m]
real, intent(in) :: temp(size(z)) ! 温度 [K]
real, intent(in) :: pres(size(z)) ! 圧力 [Pa]
real, intent(in) :: qv(size(z)) ! 水蒸気混合比 [kg/kg]
integer :: i, j, k, nz, iz_ref
real :: TLCL
real :: T_ref, P_ref, qv_ref
call interpo_search_1d( z, z_ref, iz_ref )
call interpolation_1d( z(iz_ref), z(iz_ref), temp(iz_ref), temp(iz_ref+1), z_ref, T_ref )
call interpolation_1d( z(iz_ref), z(iz_ref), pres(iz_ref), pres(iz_ref+1), z_ref, P_ref )
call interpolation_1d( z(iz_ref), z(iz_ref), qv(iz_ref), qv(iz_ref+1), z_ref, qv_ref )
TLCL=TqvP_2_TLCL( T_ref, P_ref, qv_ref ) ! z_ref のパーセルの LCL での温度.
z_LCL=z_ref+(Cpd/g)*(T_ref-TLCL)
return
end function
| Function : | |||
| z_LFC : | real | ||
| z_ref : | real, intent(in)
| ||
| z(:) : | real, intent(in)
| ||
| temp(size(z)) : | real, intent(in)
| ||
| pres(size(z)) : | real, intent(in)
| ||
| qv(size(z)) : | real, intent(in)
|
LFC 高度を計算する. ある基準高度での相当温位を上昇させたとき, 周囲の飽和相当温位と一致する高度が LFC である. ここでは, 周囲の飽和相当温位と交差する 2 点を求め, そこから高度方向に 内挿することで高度を決定する.
real function z_LFC( z_ref, z, temp, pres, qv ) ! LFC 高度を計算する.
! ある基準高度での相当温位を上昇させたとき, 周囲の飽和相当温位と一致する高度が
! LFC である.
! ここでは, 周囲の飽和相当温位と交差する 2 点を求め, そこから高度方向に
! 内挿することで高度を決定する.
use Statistics
use Thermo_Function
implicit none
real, intent(in) :: z_ref ! 基準高度 [m]
real, intent(in) :: z(:) ! 高度座標 [m]
real, intent(in) :: temp(size(z)) ! 温度 [K]
real, intent(in) :: pres(size(z)) ! 圧力 [Pa]
real, intent(in) :: qv(size(z)) ! 水蒸気混合比 [kg/kg]
integer :: i, j, k, nz, iz, iept
real :: sept(size(z))
real :: eptiz, eptiz1, ept_ref, z_sol
nz=size(z)
!-- 鉛直方向の飽和相当温位を計算.
do i=1,nz
sept(i)=thetaes_Bolton( temp(i), pres(i) )
end do
!-- z_ref での相当温位を計算する.
call interpo_search_1d( z, z_ref, iz )
!-- iz と iz+1 での相当温位を計算.
eptiz=thetae_Bolton( temp(iz), qv(iz), pres(iz) )
eptiz1=thetae_Bolton( temp(iz+1), qv(iz+1), pres(iz+1) )
!-- これらの間から相当温位を内挿
call interpolation_1d( z(iz), z(iz+1), eptiz, eptiz1, z_ref, ept_ref )
!-- z_ref より上の高度で ept_ref と飽和相当温位の関係を計算.
do i=iz+1,nz
if((sept(i)-ept_ref)*(sept(i+1)-ept_ref)<0.0)then
iept=i
exit
end if
end do
!-- パーセルの相当温位が環境場の飽和相当温位と逆転する直前の高度が iept である.
!-- ここで, iept と iept+1 の飽和相当温位から高度を内挿する.
call interpolation_1d( sept(iept), sept(iept+1), z(iept), z(iept+1), ept_ref, z_sol )
z_LFC=z_sol
return
end function
| Function : | |||
| z_LNB : | real | ||
| z_ref : | real, intent(in)
| ||
| z(:) : | real, intent(in)
| ||
| temp(size(z)) : | real, intent(in)
| ||
| pres(size(z)) : | real, intent(in)
| ||
| qv(size(z)) : | real, intent(in)
| ||
| opt : | integer, intent(in), optional
|
LNB 高度を計算する. LFC 高度以上で再び環境場の飽和相当温位と交わる高度が LNB である. ただし, 実際の測定データでは, 測器の微小変動により, LFC の直上で LNB に 到達する場合がある. そこで, opt として, LFC を求めてから LNB を計算する方法と 上空から下向きに下ろして最初に交差する高度を LNB と定義する手法の 2 パターン用意することにする. opt = 1, LFC から求める. opt = 2 直上から計算. デフォルトでは opt = 1 を計算.
real function z_LNB( z_ref, z, temp, pres, qv, opt ) ! LNB 高度を計算する.
! LFC 高度以上で再び環境場の飽和相当温位と交わる高度が LNB である.
! ただし, 実際の測定データでは, 測器の微小変動により, LFC の直上で LNB に
! 到達する場合がある.
! そこで, opt として, LFC を求めてから LNB を計算する方法と
! 上空から下向きに下ろして最初に交差する高度を LNB と定義する手法の
! 2 パターン用意することにする.
! opt = 1, LFC から求める. opt = 2 直上から計算.
! デフォルトでは opt = 1 を計算.
use Statistics
use Thermo_Function
implicit none
real, intent(in) :: z_ref ! 基準高度 [m]
real, intent(in) :: z(:) ! 高度座標 [m]
real, intent(in) :: temp(size(z)) ! 温度 [K]
real, intent(in) :: pres(size(z)) ! 圧力 [Pa]
real, intent(in) :: qv(size(z)) ! 水蒸気混合比 [kg/kg]
integer, intent(in), optional :: opt ! 計算手法のオプション
integer :: i, j, k, nz, iz, iept, counter
real :: sept(size(z))
real :: eptiz, eptiz1, ept_ref, z_sol
nz=size(z)
!-- 鉛直方向の飽和相当温位を計算.
do i=1,nz
sept(i)=thetaes_Bolton( temp(i), pres(i) )
end do
!-- z_ref での相当温位を計算する.
call interpo_search_1d( z, z_ref, iz )
!-- iz と iz+1 での相当温位を計算.
eptiz=thetae_Bolton( temp(iz), qv(iz), pres(iz) )
eptiz1=thetae_Bolton( temp(iz+1), qv(iz+1), pres(iz+1) )
!-- これらの間から相当温位を内挿
call interpolation_1d( z(iz), z(iz+1), eptiz, eptiz1, z_ref, ept_ref )
counter=0
if(present(opt))then
if(opt==2)then
do i=nz,iz+1,-1 ! 上から下に下ろす.
if((sept(i)-ept_ref)*(sept(i-1)-ept_ref)<0.0)then
iept=i-1 ! 上から下に下ろしているので, 1 要素下のデータが iept.
exit
end if
end do
else
!-- z_ref より上の高度で ept_ref と飽和相当温位の関係を計算.
do i=iz+1,nz
if((sept(i)-ept_ref)*(sept(i+1)-ept_ref)<0.0)then
counter=counter+1
if(counter==2)then
iept=i
exit
end if
end if
end do
end if
else
do i=iz+1,nz
if((sept(i)-ept_ref)*(sept(i+1)-ept_ref)<0.0)then
counter=counter+1
if(counter==2)then
iept=i
exit
end if
end if
end do
end if
!-- パーセルの相当温位が環境場の飽和相当温位と逆転する直前の高度が iept である.
!-- ここで, iept と iept+1 の飽和相当温位から高度を内挿する.
call interpolation_1d( sept(iept), sept(iept+1), z(iept), z(iept+1), ept_ref, z_sol )
z_LNB=z_sol
return
end function