Subroutine : |
|
title : | character(*), intent(in)
|
slope : | real, intent(in)
|
t_int(2) : | real, intent(in)
: | 横軸(温度)の両端値 [K] t_int(1)=tmin, t_int(2)=tmax
|
|
p_int(2) : | real, intent(in)
: | 縦軸(気圧)の両端値 [Pa] p_int(1)=pmin, p_int(2)=pmax sondeZ
が与えられる場合は, この引数は p_int(1)=zmin, p_int(2)=zmax を入れる.
|
|
sondeT : | real, intent(in), dimension(:)
|
sondeP : | real, intent(in), dimension(size(sondeT))
|
sondeQV : | real, intent(in), dimension(size(sondeT))
|
sondeZ : | real, intent(in), dimension(size(sondeT)), optional
: | ゾンデの高度 [m] sondeP : logP で縦軸をとる. このとき, pmin, pmax, sondeP
はすべて
気圧を Pa 単位でとること.
sondeZ(optional) : geometric z で縦軸をとる.
このとき, pmin, pmax は気圧ではなく,
描画したい高度 [m] の最小最大値を与える.
z=-logP の関係になるが, それは内部で処理するので, 使用者は 気にせず, min
は最低値, max は最大値を入れれば OK. デフォルトでは, 気圧で縦軸設定.
|
|
viewx_int(2) : | real, intent(in), optional
: | viewport の x 方向の両端. viewx_int(1)=viewxmin, viewx_int(2)=viewxmax
|
|
viewy_int(2) : | real, intent(in), optional
: | viewport の y 方向の両端. viewy_int(1)=viewymin, viewy_int(2)=viewymax
デフォルトは, (0.25, 0.7, 0.2, 0.8) で固定.
|
|
undef : | real, intent(in), optional
|
no_frame : | logical, intent(in), optional
: | NewFrame を呼ばない .false. = 呼ぶ. .true. = 呼ばずに NewFig. デフォルトは
.false. — 以上, 引数
|
|
skew-T ダイアグラムを作成するルーチン 横軸と縦軸のタイトルはすべて height
と temperature で固定.
subroutine SkewT( title, slope, t_int, p_int, sondeT, sondeP, sondeQV, sondeZ, viewx_int, viewy_int, undef, no_frame )
! skew-T ダイアグラムを作成するルーチン
! 横軸と縦軸のタイトルはすべて height と temperature で固定.
use Thermo_Const
implicit none
character(*), intent(in) :: title ! 表のタイトル
real, intent(in) :: slope ! skewT の傾き具合
real, intent(in) :: t_int(2) ! 横軸(温度)の両端値 [K]
! t_int(1)=tmin, t_int(2)=tmax
real, intent(in) :: p_int(2) ! 縦軸(気圧)の両端値 [Pa]
! p_int(1)=pmin, p_int(2)=pmax
! sondeZ が与えられる場合は, この引数は
! p_int(1)=zmin, p_int(2)=zmax を入れる.
real, intent(in), dimension(:) :: sondeT ! ゾンデの温度 [K]
real, intent(in), dimension(size(sondeT)) :: sondeP ! ゾンデの気圧 [Pa]
real, intent(in), dimension(size(sondeT)) :: sondeQV ! ゾンデの混合比 [kg/kg]
real, intent(in), dimension(size(sondeT)), optional :: sondeZ ! ゾンデの高度 [m]
! sondeP : logP で縦軸をとる. このとき, pmin, pmax, sondeP はすべて
! 気圧を Pa 単位でとること.
! sondeZ(optional) : geometric z で縦軸をとる.
! このとき, pmin, pmax は気圧ではなく,
! 描画したい高度 [m] の最小最大値を与える.
! z=-logP の関係になるが, それは内部で処理するので, 使用者は
! 気にせず, min は最低値, max は最大値を入れれば OK.
! デフォルトでは, 気圧で縦軸設定.
real, intent(in), optional :: viewx_int(2) ! viewport の x 方向の両端.
! viewx_int(1)=viewxmin, viewx_int(2)=viewxmax
real, intent(in), optional :: viewy_int(2) ! viewport の y 方向の両端.
! viewy_int(1)=viewymin, viewy_int(2)=viewymax
! デフォルトは, (0.25, 0.7, 0.2, 0.8) で固定.
real, intent(in), optional :: undef ! 未定義値, デフォルトでは設定しない.
logical, intent(in), optional :: no_frame ! NewFrame を呼ばない
! .false. = 呼ぶ. .true. = 呼ばずに NewFig.
! デフォルトは .false.
!-- 以上, 引数
real :: tmin ! 横軸(温度)の最低値 [K]
real :: tmax ! 横軸(温度)の最高値 [K]
real :: pmin ! 縦軸(気圧)の最低値 [Pa]
real :: pmax ! 縦軸(気圧)の最高値 [Pa]
real :: zmin ! 縦軸(高度)の最低値 [m]
real :: zmax ! 縦軸(高度)の最高値 [m]
real :: vxmin ! viewport の xmin.
real :: vxmax ! viewport の xmin.
real :: vymin ! viewport の xmin.
real :: vymax ! viewport の xmin.
!-- 以上, 引数の置き換え用変数
integer, parameter :: nx=100, ny=100
real, dimension(nx,ny) :: temp, pt, qvs, spt
real, dimension(size(sondeT)) :: sondedT, sondedP, sondedTD
real, dimension(nx) :: tempd
real, dimension(ny) :: p, z
integer :: i, j, k, nt, i_z, i_p
real :: pdmax, pdmin
real :: dp, dt
real :: undeff
real :: slope_param, dz
logical :: heiopt_z, no_frame_flag
heiopt_z=.false.
nt=size(sondeT)
if(present(undef))then
undeff=undef
else
undeff=-999.0
end if
!-- DCL 前設定
! call SGLSET('LCLIP',.true.) ! クリッピング処理
call DclSetParm( 'ENABLE_CONTOUR_MESSAGE', .false. )
!-- 以降, オプション処理
if(present(viewx_int))then
vxmin=viewx_int(1)
vxmax=viewx_int(2)
else
vxmin=0.25
vxmax=0.75
end if
if(present(viewy_int))then
vymin=viewy_int(1)
vymax=viewy_int(2)
else
vymin=0.1
vymax=0.8
end if
if(present(no_frame))then
no_frame_flag=no_frame
else
no_frame_flag=.false.
end if
!-- 引数の置き換え用変数に置き換え
pmin=p_int(1)
pmax=p_int(2)
tmin=t_int(1)
tmax=t_int(2)
dt=(tmax-tmin)/(nx-1)
tempd=(/((tmin+dt*(i-1)),i=1,nx)/)
if(present(sondeZ))then
! 高度計算した場合, z を p に直さないと pt の計算不可能 (要修正)
zmin=pmin
zmax=pmax
call interpo_search_1d( sondeZ, zmin, i_z )
call interpolation_1d( sondeZ(i_z:i_z+1), sondeP(i_z:i_z+1), zmin, pmin )
call interpo_search_1d( sondeZ, zmax, i_z )
call interpolation_1d( sondeZ(i_z:i_z+1), sondeP(i_z:i_z+1), zmax, pmax )
heiopt_z=.true.
dz=(zmax-zmin)/(ny-1)
z=(/((zmin+dz*(i-1)),i=1,ny)/)
slope_param=1.0/(zmax-zmin) ! 傾き計算の調節係数.
! だいたい, 最大高度で 1 となるように規格化
do i=1,ny
call interpo_search_1d( sondeZ, z(i), i_z )
call interpolation_1d( z(i_z:i_z+1), sondeP(i_z:i_z+1), z(i), p(i) )
end do
pdmin=zmin
pdmax=zmax
else
dp=(pmax-pmin)/(ny-1)
p=(/((pmax-dp*(i-1)),i=1,ny)/)
pdmin=pmin*1.0e-2
pdmax=pmax*1.0e-2
end if
!-- 各曲線の計算
do j=1,ny
do i=1,nx
if(heiopt_z.eqv..true.)then
temp(i,j)=tempd(i)-slope*z(j)*slope_param
else
temp(i,j)=tempd(i)+slope*log(p(j)/p0)
end if
pt(i,j)=theta_dry( temp(i,j), p(j) )
spt(i,j)=thetaes_Bolton( temp(i,j), p(j) )
qvs(i,j)=TP_2_qvs( temp(i,j), p(j) )*1.0e3
end do
end do
do j=1,nt
if(sondeT(j)==undeff.or.sondeP(j)==undeff.or.sondeQV(j)==undeff)then
sondedT(j)=undeff
sondedTD(j)=undeff
else
if(heiopt_z.eqv..true.)then
sondedT(j)=sondeT(j)+slope*sondeZ(j)*slope_param ! 2000 で割る根拠はない.
sondedTD(j)=es_TD(qvP_2_e(sondeQV(j),sondeP(j)))+slope*sondeZ(j)*slope_param
else
sondedT(j)=sondeT(j)-slope*log(sondeP(j)/p0) ! temp, tempd の操作と同じ.
sondedTD(j)=es_TD(qvP_2_e(sondeQV(j),sondeP(j)))-slope*log(sondeP(j)/p0)
end if
end if
if(heiopt_z.eqv..true.)then
sondedP(j)=sondeZ(j)
else
sondedP(j)=sondeP(j)*1.0e-2
end if
end do
if(no_frame_flag.eqv..true.)then
call DclNewFig
else
call DclNewFrame
end if
call DclScalingPoint( sondedT, sondedP )
if(heiopt_z.eqv..true.)then
call DclSetTransNumber( 1 )
call DclSetWindow( tmin, tmax, pdmin, pdmax )
else
call DclSetTransNumber( 2 )
call DclSetWindow( tmin, tmax, pdmax, pdmin )
end if
call DclSetViewPort( vxmin, vxmax, vymin, vymax )
call DclSetTransFunction
if(heiopt_z.eqv..true.)then
call DclSetTitle( 'Temperature', 'height', '', '' )
else
call DclSetTitle( 'Temperature', 'pressure', '', '' )
end if
call DclDrawScaledAxis
! call DclDrawAxis( 'v', pdmax/40.0, pdmax/10.0, section=vxmin )
! call DclDrawAxis( 'v', pdmax/40.0, pdmax/10.0, section=vxmax )
call UDISET('INDXMJ', 31)
call UDISET('INDXMN', 31)
call DclDrawContour( temp )
call UDISET('INDXMJ', 22)
call UDISET('INDXMN', 22)
call DclDrawContour( pt )
call UDISET('INDXMJ', 13)
call UDISET('INDXMN', 13)
call DclDrawContour( spt )
call UDISET('INDXMJ', 1)
call DclDrawContour( qvs )
call DclDrawLine( sondedT, sondedP, index=43 )
call DclDrawLine( sondedTD, sondedP, index=43, type=2 )
call DclDrawTitle( 't', trim(title) )
end subroutine