Class Dcl_Thermo_Graph
In: dcl_thermo_graph.f90

熱力学関係のグラフの自動化ルーチン

Methods

SkewT   SkewT_vec  

Included Modules

Thermo_Const Thermo_Function Statistics dcl

Public Instance methods

Subroutine :
title :character(*), intent(in)
: 表のタイトル
slope :real, intent(in)
: skewT の傾き具合
tmin :real, intent(in)
: 横軸(温度)の最低値 [K]
tmax :real, intent(in)
: 横軸(温度)の最高値 [K]
pmin :real, intent(in)
: 縦軸(気圧)の最低値 [Pa]
pmax :real, intent(in)
: 縦軸(気圧)の最高値 [Pa]
sondeT :real, intent(in), dimension(:)
: ゾンデの温度 [K]
sondeP :real, intent(in), dimension(size(sondeT))
: ゾンデの気圧 [Pa]
sondeQV :real, intent(in), dimension(size(sondeT))
: ゾンデの混合比 [kg/kg]
heiopt :character(1), intent(in), optional
: 縦軸のとり方 heiopt=‘p’ : logP で縦軸をとる. このとき, pmin, pmax, sondeP はすべて
             気圧を Pa 単位でとること.

heiopt=‘g’ : geometric z で縦軸をとる. このとき, pmin, pmax, sondeP は

             すべて m で単位をとること.

z=-logP の関係になるが, それは内部で処理するので, 使用者は 気にせず, min は最低値, max は最大値を入れれば OK. デフォルトでは, 気圧で縦軸設定.

viewxmin :real, intent(in), optional
: viewport の xmin.
viewxmax :real, intent(in), optional
: viewport の xmin.
viewymin :real, intent(in), optional
: viewport の xmin.
viewymax :real, intent(in), optional
: viewport の xmin. デフォルトは, (0.25, 0.7, 0.2, 0.8) で固定.
undef :real, intent(in), optional
: 未定義値, デフォルトでは設定しない.

skew-T ダイアグラムを作成するルーチン 横軸と縦軸のタイトルはすべて height と temperature で固定.

[Source]

subroutine SkewT( title, slope, tmin, tmax, pmin, pmax, sondeT, sondeP, sondeQV, heiopt, viewxmin, viewxmax, viewymin, viewymax, undef )
! skew-T ダイアグラムを作成するルーチン
! 横軸と縦軸のタイトルはすべて height と temperature で固定.
  use Thermo_Const
  use Thermo_Function
  use Dcl
  implicit none
  character(*), intent(in) :: title  ! 表のタイトル
  real, intent(in) :: slope  ! skewT の傾き具合
  real, intent(in) :: tmin  ! 横軸(温度)の最低値 [K]
  real, intent(in) :: tmax  ! 横軸(温度)の最高値 [K]
  real, intent(in) :: pmin  ! 縦軸(気圧)の最低値 [Pa]
  real, intent(in) :: pmax  ! 縦軸(気圧)の最高値 [Pa]
  real, intent(in), dimension(:) :: sondeT  ! ゾンデの温度 [K]
  real, intent(in), dimension(size(sondeT)) :: sondeP  ! ゾンデの気圧 [Pa]
  real, intent(in), dimension(size(sondeT)) :: sondeQV  ! ゾンデの混合比 [kg/kg]
  character(1), intent(in), optional :: heiopt  ! 縦軸のとり方
  ! heiopt='p' : logP で縦軸をとる. このとき, pmin, pmax, sondeP はすべて
  !              気圧を Pa 単位でとること.
  ! heiopt='g' : geometric z で縦軸をとる. このとき, pmin, pmax, sondeP は
  !              すべて m で単位をとること.
  ! z=-logP の関係になるが, それは内部で処理するので, 使用者は
  ! 気にせず, min は最低値, max は最大値を入れれば OK.
  ! デフォルトでは, 気圧で縦軸設定.
  real, intent(in), optional :: viewxmin  ! viewport の xmin.
  real, intent(in), optional :: viewxmax  ! viewport の xmin.
  real, intent(in), optional :: viewymin  ! viewport の xmin.
  real, intent(in), optional :: viewymax  ! viewport の xmin.
  ! デフォルトは, (0.25, 0.7, 0.2, 0.8) で固定.
  real, intent(in), optional :: undef  ! 未定義値, デフォルトでは設定しない.
  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
  integer :: i, j, k, nt
  real :: pdmax, pdmin
  real :: dp, dpd, dt
  real :: vxmin, vxmax, vymin, vymax, undeff

  nt=size(sondeT)

  pdmin=pmin*1.0e-2
  pdmax=pmax*1.0e-2

  if(present(undef))then
     undeff=undef
  end if

  if(present(viewxmin))then
     vxmin=viewxmin
  else
     vxmin=0.25
  end if
  if(present(viewxmax))then
     vxmax=viewxmax
  else
     vxmax=0.75
  end if
  if(present(viewymin))then
     vymin=viewymin
  else
     vymin=0.1
  end if
  if(present(viewymax))then
     vymax=viewymax
  else
     vymax=0.8
  end if

  dt=(tmax-tmin)/(nx-1)
  tempd=(/((tmin+dt*(i-1)),i=1,nx)/)

  if(present(heiopt))then
     if(heiopt=='p')then
        dp=(pmax-pmin)/(ny-1)
        dpd=(pdmax-pdmin)/(ny-1)
        p=(/((pmax-dp*(i-1)),i=1,ny)/)
     else  ! 高度計算した場合, z を p に直さないと pt の計算不可能 (要修正)
        dp=(pmax-pmin)/(ny-1)
        dpd=(pdmax-pdmin)/(ny-1)
        p=(/((pmin+dp*(i-1)),i=1,ny)/)
     end if
  else
     dp=(pmax-pmin)/(ny-1)
     dpd=(pdmax-pdmin)/(ny-1)
     p=(/((pmax-dp*(i-1)),i=1,ny)/)
  end if

!-- 各曲線の計算
  do j=1,ny
     do i=1,nx
        temp(i,j)=tempd(i)+slope*log(p(j)/p0)
        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
        sondedP(j)=undeff
        sondedTD(j)=undeff
     else
        sondedT(j)=sondeT(j)-slope*log(sondeP(j)/p0)  ! temp, tempd の操作と同じ.
        sondedP(j)=sondeP(j)*1.0e-2
        sondedTD(j)=es_TD(qvP_2_e(sondeQV(j),sondeP(j)))-slope*log(sondeP(j)/p0)
     end if
  end do

  call DclNewFrame

  call DclScalingPoint( sondedT, sondedP )
  call DclSetTransNumber( 2 )
  call DclSetWindow( tmin, tmax, pdmax, pdmin )
  call DclSetViewPort( vxmin, vxmax, vymin, vymax )
  call DclSetTransFunction

  call DclSetTitle( 'Temperature', 'height', '', '' )

  call DclDrawScaledAxis

  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 :
title :character(*), intent(in)
: 表のタイトル
slope :real, intent(in)
: skewT の傾き具合
tmin :real, intent(in)
: 横軸(温度)の最低値 [K]
tmax :real, intent(in)
: 横軸(温度)の最高値 [K]
pmin :real, intent(in)
: 縦軸(気圧)の最低値 [Pa]
pmax :real, intent(in)
: 縦軸(気圧)の最高値 [Pa]
sondeT :real, intent(in), dimension(:)
: ゾンデの温度 [K]
sondeP :real, intent(in), dimension(size(sondeT))
: ゾンデの気圧 [Pa]
sondeQV :real, intent(in), dimension(size(sondeT))
: ゾンデの混合比 [kg/kg]
sondeU :real, intent(in), dimension(size(sondeT))
: ゾンデの東西風 [m/s]
sondeV :real, intent(in), dimension(size(sondeT))
: ゾンデの南北風 [m/s]
heiopt :character(1), intent(in), optional
: 縦軸のとり方 heiopt=‘p’ : logP で縦軸をとる. このとき, pmin, pmax, sondeP はすべて
             気圧を Pa 単位でとること.

heiopt=‘g’ : geometric z で縦軸をとる. このとき, pmin, pmax, sondeP は

             すべて m で単位をとること.

z=-logP の関係になるが, それは内部で処理するので, 使用者は 気にせず, min は最低値, max は最大値を入れれば OK. デフォルトでは, 気圧で縦軸設定.

viewxmin :real, intent(in), optional
: viewport の xmin.
viewxmax :real, intent(in), optional
: viewport の xmin.
viewymin :real, intent(in), optional
: viewport の xmin.
viewymax :real, intent(in), optional
: viewport の xmin. デフォルトは, (0.25, 0.7, 0.2, 0.8) で固定.
vx :integer, intent(in), optional
: ベクトルの間引き格子数. 間引いた場合, 圧力に関して等格子上に配置する. デフォルトは 10 格子
undef :real, intent(in), optional
: 未定義値, デフォルトでは設定しない.

skew-T ダイアグラム + vector 図を作成するルーチン 横軸と縦軸のタイトルはすべて height と temperature で固定.

[Source]

subroutine SkewT_vec( title, slope, tmin, tmax, pmin, pmax, sondeT, sondeP, sondeQV, sondeU, sondeV, heiopt, viewxmin, viewxmax, viewymin, viewymax, vx, undef )
! skew-T ダイアグラム + vector 図を作成するルーチン
! 横軸と縦軸のタイトルはすべて height と temperature で固定.
  use Thermo_Const
  use Thermo_Function
  use Statistics
  use Dcl
  implicit none
  character(*), intent(in) :: title  ! 表のタイトル
  real, intent(in) :: slope  ! skewT の傾き具合
  real, intent(in) :: tmin  ! 横軸(温度)の最低値 [K]
  real, intent(in) :: tmax  ! 横軸(温度)の最高値 [K]
  real, intent(in) :: pmin  ! 縦軸(気圧)の最低値 [Pa]
  real, intent(in) :: pmax  ! 縦軸(気圧)の最高値 [Pa]
  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)) :: sondeU  ! ゾンデの東西風 [m/s]
  real, intent(in), dimension(size(sondeT)) :: sondeV  ! ゾンデの南北風 [m/s]
  character(1), intent(in), optional :: heiopt  ! 縦軸のとり方
  ! heiopt='p' : logP で縦軸をとる. このとき, pmin, pmax, sondeP はすべて
  !              気圧を Pa 単位でとること.
  ! heiopt='g' : geometric z で縦軸をとる. このとき, pmin, pmax, sondeP は
  !              すべて m で単位をとること.
  ! z=-logP の関係になるが, それは内部で処理するので, 使用者は
  ! 気にせず, min は最低値, max は最大値を入れれば OK.
  ! デフォルトでは, 気圧で縦軸設定.
  real, intent(in), optional :: viewxmin  ! viewport の xmin.
  real, intent(in), optional :: viewxmax  ! viewport の xmin.
  real, intent(in), optional :: viewymin  ! viewport の xmin.
  real, intent(in), optional :: viewymax  ! viewport の xmin.
  ! デフォルトは, (0.25, 0.7, 0.2, 0.8) で固定.
  integer, intent(in), optional :: vx  ! ベクトルの間引き格子数.
  ! 間引いた場合, 圧力に関して等格子上に配置する.
  ! デフォルトは 10 格子
  real, intent(in), optional :: undef  ! 未定義値, デフォルトでは設定しない.
  integer, parameter :: nx=100, ny=100
  real, dimension(nx,ny) :: temp, pt, qvs, spt
  real, dimension(size(sondeT)) :: sondedT, sondedP, sondedTD, sondelogP
  real, dimension(nx) :: tempd
  real, dimension(ny) :: p
  integer :: i, j, k, nt, vecx
  real :: pdmax, pdmin
  real :: dp, dpd, dt, dpv
  real :: vxmin, vxmax, vymin, vymax, undeff
  real, allocatable, dimension(:,:) :: u, v
  real, allocatable, dimension(:) :: pv
  integer, allocatable, dimension(:) :: interpo

  nt=size(sondeT)

  pdmin=pmin*1.0e-2
  pdmax=pmax*1.0e-2

  if(present(vx))then
     vecx=vx
  else
     vecx=10
  end if

  if(present(undef))then
     undeff=undef
  end if

  if(present(viewxmin))then
     vxmin=viewxmin
  else
     vxmin=0.25
  end if
  if(present(viewxmax))then
     vxmax=viewxmax
  else
     vxmax=0.75
  end if
  if(present(viewymin))then
     vymin=viewymin
  else
     vymin=0.1
  end if
  if(present(viewymax))then
     vymax=viewymax
  else
     vymax=0.8
  end if

  dt=(tmax-tmin)/(nx-1)
  tempd=(/((tmin+dt*(i-1)),i=1,nx)/)

  if(present(heiopt))then
     if(heiopt=='p')then
        dp=(pmax-pmin)/(ny-1)
        dpd=(pdmax-pdmin)/(ny-1)
        p=(/((pmax-dp*(i-1)),i=1,ny)/)
     else  ! 高度計算した場合, z を p に直さないと pt の計算不可能 (要修正)
        dp=(pmax-pmin)/(ny-1)
        dpd=(pdmax-pdmin)/(ny-1)
        p=(/((pmin+dp*(i-1)),i=1,ny)/)
     end if
  else
     dp=(pmax-pmin)/(ny-1)
     dpd=(pdmax-pdmin)/(ny-1)
     p=(/((pmax-dp*(i-1)),i=1,ny)/)
  end if

!-- 各曲線の計算
  do j=1,ny
     do i=1,nx
        temp(i,j)=tempd(i)+slope*log(p(j)/p0)
        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.or.sondeU(j)==undeff.or.sondeV(j)==undeff)then
        sondedT(j)=undeff
        sondedP(j)=undeff
        sondedTD(j)=undeff
        sondelogP(j)=undeff
     else
        sondedT(j)=sondeT(j)-slope*log(sondeP(j)/p0)  ! temp, tempd の操作と同じ.
        sondedP(j)=sondeP(j)*1.0e-2
        sondedTD(j)=es_TD(qvP_2_e(sondeQV(j),sondeP(j)))-slope*log(sondeP(j)/p0)
        sondelogP(j)=-log(sondeP(j))
     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
     pv(i)=-log(pmax-dpv*(i-1))  ! 高度として代入.
     if(sondelogP(nt)>pv(i))then
        call interpo_search_1d( sondelogP, pv(i), interpo(i), undeff )
        
        if(sondeU(interpo(i))/=undeff.and.sondeU(interpo(i)+1)/=undeff.and.sondeV(interpo(i))/=undeff.and.sondeV(interpo(i)+1)/=undeff)then
!        call interpolation_1d( sondeZ(interpo(i)), sondeZ(interpo(i)+1),  &
!  &          sondeU(interpo(i)), sondeU(interpo(i)+1), pv(i), u(2,i))
!        call interpolation_1d( sondelogP(interpo(i)), sondelogP(interpo(i)+1),  &
!  &          sondeV(interpo(i)), sondeV(interpo(i)+1), pv(i), v(2,i))
!write(*,*) sondeU(interpo(i)), sondeU(interpo(i)+1)
! 本当は上のコメントアウトのように内挿で正確に計算すべきであるが,
! データ点が細かいことから, もっとも近い点で近似してもよいとする.
           u(2,i)=sondeU(interpo(i))
           v(2,i)=sondeV(interpo(i))
        else
           u(2,i)=undeff
           v(2,i)=undeff
        end if
     else
        u(2,i)=undeff
        v(2,i)=undeff
     end if
  end do

  call DclNewFrame

  call DclScalingPoint( sondedT, sondedP )
  call DclSetTransNumber( 2 )
  call DclSetWindow( tmin, tmax, pdmax, pdmin )
  call DclSetViewPort( vxmin, vxmax, vymin, vymax )
  call DclSetTransFunction

  call DclSetTitle( 'Temperature', 'height', '', '' )

  call DclDrawScaledAxis

  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) )

  call GRFIG
  call DclSetWindow( 0.0, 1.0, pdmax, pdmin )
  call DclSetViewPort( vxmax, vxmax+0.1, vymin, vymax )
  call DclSetTransNumber( 2 )
  call DclSetTransFunction

  CALL UZLSET( 'LABELYR', .FALSE. )
  CALL UZLSET( 'LABELYL', .FALSE. )

  call DclDrawAxis( 'v', pdmax/40.0, pdmax/10.0, section=vxmax+0.1 )

  call UGLSET( 'LMSG', .false. )

  call DclDrawVectors( u, v )
write(*,*) u(2,:)

end subroutine