| Class | Dcl_Thermo_Graph |
| In: |
dcl_thermo_graph.f90
|
熱力学関係のグラフの自動化ルーチン
| Subroutine : | |||
| title : | character(*), intent(in)
| ||
| slope : | real, intent(in)
| ||
| t_int(2) : | real, intent(in)
| ||
| p_int(2) : | real, intent(in)
| ||
| 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
| ||
| viewx_int(2) : | real, intent(in), optional
| ||
| viewy_int(2) : | real, intent(in), optional
| ||
| undef : | real, intent(in), optional
| ||
| no_frame : | logical, intent(in), optional
|
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
| Subroutine : | |||
| p_int(2) : | real, intent(in)
| ||
| vx : | real, intent(in)
| ||
| sondeP : | real, intent(in), dimension(:)
| ||
| sondeU : | real, intent(in), dimension(size(sondeP))
| ||
| sondeV : | real, intent(in), dimension(size(sondeP))
| ||
| sondeZ : | real, intent(in), dimension(size(sondeP)), optional
| ||
| view_int(2) : | real, intent(in), optional
| ||
| vn : | integer, intent(in), optional
| ||
| undef : | real, intent(in), optional
| ||
| unit_flag : | logical, intent(in), optional
| ||
| unity : | real, intent(in), optional
| ||
| u_posi(2) : | real, intent(in), optional
| ||
| uni_len : | real, intent(in), optional
|
skew-T ダイアグラム + vector 図を作成するルーチン 横軸と縦軸のタイトルはすべて height と temperature で固定.
subroutine axis_vec( p_int, vx, sondeP, sondeU, sondeV, sondeZ, view_int, vn, undef, unit_flag, unity, u_posi, uni_len )
! skew-T ダイアグラム + vector 図を作成するルーチン
! 横軸と縦軸のタイトルはすべて height と temperature で固定.
use Thermo_Const
use Thermo_Function
use Statistics
use Dcl
implicit none
real, intent(in) :: p_int(2) ! 横軸(気圧)の両端値 [Pa]
! p_int(1)=pmin, p_int(2)=pmax
real, intent(in) :: vx ! 軸線を描く vx 座標.
real, intent(in), dimension(:) :: sondeP ! ゾンデの気圧 [Pa]
real, intent(in), dimension(size(sondeP)) :: sondeU ! ゾンデの東西風 [m/s]
real, intent(in), dimension(size(sondeP)) :: sondeV ! ゾンデの南北風 [m/s]
real, intent(in), dimension(size(sondeP)), 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 :: view_int(2) ! viewport の両端.
! view_int(1)=vmin, view_int(2)=vmax
! デフォルトは, (0.2, 0.8) で固定.
integer, intent(in), optional :: vn ! ベクトルの間引き格子数.
! 間引いた場合, 描画座標に関して等格子上に配置する.
! デフォルトは 10 格子
real, intent(in), optional :: undef ! 未定義値, デフォルトでは設定しない.
logical, intent(in), optional :: unit_flag ! 単位ベクトルを描くか.
! default = .false. : 描かない.
real, intent(in), optional :: unity ! unit_flag = .true. の時,
! 単位ベクトルの大きさ.
! 指定しない場合, エラーで終了.
! V 座標系では 0.1 で一定.
real, intent(in), optional :: u_posi(2) ! 単位ベクトルの原点座標 (V 座標).
! u_posi=(/vx,vy/)
! 指定しない場合, 下端の -0.1 のところに描く.
real, intent(in), optional :: uni_len ! 単位ベクトルの長さ (V 座標).
! default = 0.1
!-- 以上, 引数
real :: pmin ! 縦軸(気圧)の最低値 [Pa]
real :: pmax ! 縦軸(気圧)の最高値 [Pa]
real :: vmin ! viewport の min.
real :: vmax ! viewport の max.
logical :: u_flag
!-- 以上, 引数の置き換え用変数
real, dimension(size(sondeP)) :: sondedP, sondelogP
integer :: i, j, k, nt, vecx
real :: pdmax, pdmin, dpv
real :: undeff
real :: uvpx, uvpy, ulen
real, allocatable, dimension(:,:) :: u, v
real, allocatable, dimension(:) :: pv
integer, allocatable, dimension(:) :: interpo
logical :: heiopt_z, no_frame_flag
heiopt_z=.false.
nt=size(sondeP)
if(present(undef))then
undeff=undef
else
undeff=-999.0
end if
if(present(vn))then
vecx=vn
else
vecx=10
end if
!-- 以降, オプション処理
if(present(view_int))then
vmin=view_int(1)
vmax=view_int(2)
else
vmin=0.2
vmax=0.8
end if
if(present(unit_flag))then
u_flag=unit_flag
if(u_flag.eqv..true.)then
if(present(unity))then
write(*,*) "unity = ", unity
else
write(*,*) "*** ERROR (axis_vec) ***"
write(*,*) "unit_flag is true, so unity must be set."
write(*,*) "STOP."
stop
end if
end if
else
u_flag=.false.
end if
if(present(u_posi))then
uvpx=u_posi(1)
uvpy=u_posi(2)
else
uvpx=vx
uvpy=vmin-0.1
end if
if(present(uni_len))then
ulen=uni_len
else
ulen=0.1
end if
!-- 引数の置き換え用変数に置き換え
pmin=p_int(1)
pmax=p_int(2)
if(present(sondeZ))then
! 高度計算した場合, z を p に直さないと pt の計算不可能 (要修正)
pdmax=pmax
pdmin=pmin
heiopt_z=.true.
else
pdmin=pmin*1.0e-2
pdmax=pmax*1.0e-2
end if
do j=1,nt
if(present(sondeZ))then
sondelogP(j)=sondeZ(j) ! おそらく, この値を undef にすると DCL でエラーとなる
else
sondelogP(j)=-log(sondeP(j)) ! おそらく, この値を undef にすると DCL でエラーとなる
end if
end do
!-- ベクトルの間引き
allocate(pv(vecx))
allocate(interpo(vecx))
allocate(u(2,vecx))
allocate(v(2,vecx))
u(1,1:vecx)=0.0
v(1,1:vecx)=0.0
dpv=(pmax-pmin)/(vecx-1)
do i=1,vecx
u(2,i)=undeff
v(2,i)=undeff
if(present(sondeZ))then
pv(i)=pmin+dpv*(i-1) ! 高度として代入.
else
pv(i)=-log(pmax-dpv*(i-1)) ! 高度として代入.
end if
if(sondelogP(nt)>pv(i))then
call interpo_search_1d( sondelogP, pv(i), interpo(i), int(undeff) )
if(interpo(i)/=int(undeff))then
if(sondeU(interpo(i))/=undeff.and.sondeU(interpo(i)+1)/=undeff.and. sondeV(interpo(i))/=undeff.and.sondeV(interpo(i)+1)/=undeff)then
! 本当は上のコメントアウトのように内挿で正確に計算すべきであるが,
! データ点が細かいことから, もっとも近い点で近似してもよいとする.
u(2,i)=sondeU(interpo(i))
v(2,i)=sondeV(interpo(i))
end if
end if
end if
end do
!-- DCL 前設定
! call SGLSET('LCLIP',.false.) ! クリッピング処理
call DclNewFig
if(heiopt_z.eqv..true.)then
call DclSetTransNumber( 1 )
call DclSetWindow( 0.0, 1.0, pdmin, pdmax )
else
call DclSetTransNumber( 2 )
call DclSetWindow( 0.0, 1.0, pdmax, pdmin )
end if
call DclSetViewPort( vx, vx+0.1, vmin, vmax )
call DclSetTransFunction
CALL UZLSET( 'LABELYR', .FALSE. )
CALL UZLSET( 'LABELYL', .FALSE. )
call DclDrawAxis( 'v', pdmax/40.0, pdmax/10.0, section=vx )
call UGLSET( 'LMSG', .false. )
call DclDrawVectors( u, v )
!-- ベクトルに関する処理
if(u_flag.eqv..true.)then
call DclSetParm( 'VECTOR:LUMSG', .false. )
!-- ベクトルの V 座標系変換処理
call DclSetParm( 'VECTOR:LNRMAL', .false. )
call DclSetParm( 'VECTOR:XFACT1', ulen/unity )
call DclSetParm( 'VECTOR:YFACT1', ulen/unity )
!-- 以下, 単位ベクトル用ルーチンで描こうとしたが, うまく描けない.
!-- 単位ベクトルを描くか.
! call DclSetParm( 'VECTOR:LUNIT', .true. )
!-- 単位ベクトルの長さ, 描画原点指定.
! call DclSetParm( 'VECTOR:VXUNIT', ulen )
! call DclSetParm( 'VECTOR:VYUNIT', 0.0 )
! call DclSetParm( 'VECTOR:VXULOC', uvpx )
! call DclSetParm( 'VECTOR:VYULOC', uvpy )
! call DclSetUnitVectorTitle( 'X', 'U='// &
! & trim(adjustl(r2c_convert(unity,'(f6.2)')))//'m s|-1"' )
!-- ので, ad hoc.
call DclDrawArrowNormalized( uvpx, uvpy, uvpx+ulen, uvpy, index=3 )
call DclDrawTextNormalized( uvpx+0.5*ulen, uvpy-0.025, 'U='// trim(adjustl(r2c_convert(unity,'(f6.2)')))//'m s|-1"', height=0.015, index=3 )
end if
end subroutine axis_vec