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 の傾き具合
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
sondeT :real, intent(in), dimension(:)
: ゾンデの温度 [K]
sondeP :real, intent(in), dimension(size(sondeT))
: ゾンデの気圧 [Pa]
sondeQV :real, intent(in), dimension(size(sondeT))
: ゾンデの混合比 [kg/kg]
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 で固定.

[Source]

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
  use Thermo_Function
  use Statistics
  use Dcl
  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
  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, dpd, 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
  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
  else
     dp=(pmax-pmin)/(ny-1)
     p=(/((pmax-dp*(i-1)),i=1,ny)/)
  end if

  pdmin=pmin*1.0e-2
  pdmax=pmax*1.0e-2
  dpd=(pdmax-pdmin)/(ny-1)

!-- 各曲線の計算
  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, zmin, zmax )
  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 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 の傾き具合
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
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]
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) で固定.
vx :integer, intent(in), optional
: ベクトルの間引き格子数. 間引いた場合, 圧力に関して等格子上に配置する. デフォルトは 10 格子
undef :real, intent(in), optional
: 未定義値, デフォルトでは設定しない.
no_frame :logical, intent(in), optional
: NewFrame を呼ばない .false. = 呼ぶ. .true. = 呼ばずに NewFig. デフォルトは .false. — 以上, 引数

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

[Source]

subroutine SkewT_vec( title, slope, t_int, p_int, sondeT, sondeP, sondeQV, sondeU, sondeV, sondeZ, viewx_int, viewy_int, vx, undef, no_frame )
! 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) :: 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
  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]
  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) で固定.
  integer, intent(in), optional :: vx  ! ベクトルの間引き格子数.
  ! 間引いた場合, 圧力に関して等格子上に配置する.
  ! デフォルトは 10 格子
  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 :: 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, sondelogP
  real, dimension(nx) :: tempd
  real, dimension(ny) :: p
  integer :: i, j, k, nt, vecx
  real :: pdmax, pdmin
  real :: dp, dpd, dt, dpv
  real :: undeff
  real, allocatable, dimension(:,:) :: u, v
  real, allocatable, dimension(:) :: pv
  integer, allocatable, dimension(:) :: interpo
  logical :: heiopt_z, no_frame_flag

  heiopt_z=.false.
  nt=size(sondeT)

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

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

!-- DCL 前設定
  call SGLSET('LCLIP',.true.)  ! クリッピング処理

!-- 以降, オプション処理
  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)

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

!-- DCL 前設定
  call SGLSET('LCLIP',.true.)  ! クリッピング処理
  call DclSetParm( 'ENABLE_CONTOUR_MESSAGE', .false. )

  dt=(tmax-tmin)/(nx-1)
  tempd=(/((tmin+dt*(i-1)),i=1,nx)/)
  dp=(pmax-pmin)/(ny-1)
  p=(/((pmax-dp*(i-1)),i=1,ny)/)

  if(present(sondeZ))then
     ! 高度計算した場合, z を p に直さないと pt の計算不可能 (要修正)
     if(present(sondeZ))then
        dpd=(pmax-pmin)/(ny-1)
        pdmax=pmax
        pdmin=pmin
        heiopt_z=.true.
     end if
  else
     pdmin=pmin*1.0e-2
     pdmax=pmax*1.0e-2
     dpd=(pdmax-pdmin)/(ny-1)
  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)=sondeP(j)*1.0e-2  ! おそらく, この値を undef にすると DCL でエラーとなる
        sondedTD(j)=undeff
        sondelogP(j)=-log(sondeP(j))  ! おそらく, この値を undef にすると DCL でエラーとなる
     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
     if(heiopt_z.eqv..true.)then
        sondedP(j)=sondeZ(j)
     else
        sondedP(j)=sondeP(j)*1.0e-2
     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), 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
!        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
     else
        u(2,i)=undeff
        v(2,i)=undeff
     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 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) )

!-- DCL 前設定
  call SGLSET('LCLIP',.false.)  ! クリッピング処理

  call GRFIG
  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( vxmax, vxmax+0.1, vymin, vymax )
  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