Class typhoon_analy
In: typhoon_analy.f90

台風感度実験用スペシャル解析モジュール

Methods

Included Modules

analy max_min statistics Geometry Thermo_Const Phys_Const Algebra

Public Instance methods

Subroutine :
r(:) :real, intent(in)
: r 方向の位置座標 [m]
coril(size(r)) :real, intent(in)
: コリオリパラメータ [/s]
v(size(r)) :real, intent(in)
: r 方向の位置座標 [m]
rho(size(r)) :real, intent(in)
: 密度 [kg/m^3]
r_ref :real, intent(in)
: 積分定数となる位置座標 [m]
p_ref :real, intent(in)
: r_ref での気圧 (積分定数) [Pa]
pres(size(r)) :real, intent(inout)
: 傾度風平衡での気圧 [kg/m^3]
 傾度風平衡場を満たす気圧場を計算する.

[Source]

subroutine grad_wind_pres( r, coril, v, rho, r_ref, p_ref, pres )
!  傾度風平衡場を満たす気圧場を計算する.
  use Algebra
  implicit none
  real, intent(in) :: r(:)  ! r 方向の位置座標 [m]
  real, intent(in) :: coril(size(r))  ! コリオリパラメータ [/s]
  real, intent(in) :: v(size(r))  ! r 方向の位置座標 [m]
  real, intent(in) :: rho(size(r))  ! 密度 [kg/m^3]
  real, intent(in) :: r_ref  ! 積分定数となる位置座標 [m]
  real, intent(in) :: p_ref  ! r_ref での気圧 (積分定数) [Pa]
  real, intent(inout) :: pres(size(r))  ! 傾度風平衡での気圧 [kg/m^3]
  integer :: i, nr
  real :: grad(size(r))

  nr=size(r)

  do i=1,nr
     if(r(i)/=0.0)then
        grad(i)=rho(i)*(v(i)*v(i)/r(i)+coril(i)*v(i))
     else
        grad(i)=0.0
     end if
  end do

  do i=1,nr
     if(r(i)<r_ref)then
        call rectangle_int( r, grad, r(i), r_ref, pres(i) )
        pres(i)=p_ref-pres(i)
     else
        if(r(i)>r_ref)then
           call rectangle_int( r, grad, r_ref, r(i), pres(i) )
           pres(i)=p_ref+pres(i)
        else
           pres(i)=p_ref
        end if
     end if
  end do

end subroutine
Subroutine :
r(:) :real, intent(in)
: 動径座標 [m]
z(:) :real, intent(in)
: 鉛直座標 [m]
coril(size(r),size(z)) :real, intent(in)
: コリオリパラメータ [/s]
v(size(r),size(z)) :real, intent(in)
: 軸対称流 [m/s]
pres_s(size(z)) :real, intent(in)
: サウンディングの気圧 [Pa]
rho_s(size(z)) :real, intent(in)
: サウンディングの密度 [kg/m^3]
pres(size(r),size(z)) :real, intent(inout)
: 平衡場の気圧 [Pa]
rho(size(r),size(z)) :real, intent(inout)
: 平衡場の密度 [K]
error :real, intent(in), optional
: イタレーションの収束条件 default = 1.0e-5
 サウンディングと軸対称流から静力学・傾度風平衡場の計算.

[Source]

subroutine hydro_grad_eqb( r, z, coril, v, pres_s, rho_s, pres, rho, error )
!  サウンディングと軸対称流から静力学・傾度風平衡場の計算.
  use Thermo_Const
  use Phys_Const
  use Algebra
  use Analy
  implicit none
  real, intent(in) :: r(:)  ! 動径座標 [m]
  real, intent(in) :: z(:)  ! 鉛直座標 [m]
  real, intent(in) :: coril(size(r),size(z))  ! コリオリパラメータ [/s]
  real, intent(in) :: v(size(r),size(z))  ! 軸対称流 [m/s]
  real, intent(in) :: pres_s(size(z))  ! サウンディングの気圧 [Pa]
  real, intent(in) :: rho_s(size(z))  ! サウンディングの密度 [kg/m^3]
  real, intent(in), optional :: error  ! イタレーションの収束条件
                    ! default = 1.0e-5
  real, intent(inout) :: pres(size(r),size(z))  ! 平衡場の気圧 [Pa]
  real, intent(inout) :: rho(size(r),size(z))  ! 平衡場の密度 [K]
  real :: old_pres(size(r),size(z)), old_rho(size(r),size(z))
  integer :: nr, nz
  integer :: i, j, k
  real :: err, err_tmp, err_max

  nr=size(r)
  nz=size(z)

  if(present(error))then
     err_max=error
  else
     err_max=1.0e-5
  end if

!-- 以下で各高度において, 密度は一定であるとして傾度風平衡から気圧を計算,
!-- その値を用いて静力学平衡から密度を修正. eps 以下になるまで繰り返す.
!-- 外縁で 2 次元場とサウンディングを一致.
  do i=1,nz
     old_pres(nr,i)=pres_s(i)
  end do
!-- 密度については, 水平面一様で設定
  do j=1,nz
     do i=1,nr
        old_rho(i,j)=rho_s(j)
     end do
  end do

!-- 以下でイタレーション開始.
  err=err_max

  do while(err>=err_max)
     err=0.0
!-- 傾度風平衡から圧力場を計算
     do j=1,nz
        call grad_wind_pres( r, coril(:,j), v(:,j), old_rho(:,j), r(nr), old_pres(nr,j), pres(:,j) )
     end do

!-- 静力学平衡から密度場を修正
     do i=1,nr
        call grad_1d( z, pres(i,:), rho(i,:) )
        do j=1,nz
if(i==1)then
write(*,*) "#### pres", pres(i,j), rho(i,j)
end if
           rho(i,j)=-rho(i,j)/g  ! 静力学の式から, dp/dz=-g*rho であるので
        end do
     end do

!-- 密度場の収束を計算
     do j=1,nz
        do i=1,nr
           if(rho(i,j)==0.0)then
              err_tmp=abs(old_rho(i,j)-rho(i,j))/abs(old_rho(i,j))
           else
              err_tmp=abs(old_rho(i,j)-rho(i,j))/abs(old_rho(i,j))
           end if

!-- 最大誤差の更新
           if(err<=err_tmp)then
              err=err_tmp
           end if

           old_rho(i,j)=rho(i,j)
           old_pres(i,j)=pres(i,j)

        end do
     end do

  end do

end subroutine
Subroutine :
signal :character(2)
: 計算するテンソル成分.
‘11’, ‘22’, ‘33‘
= それぞれ対角テンソル成分
‘12’, ‘13’, ‘21’, ‘23’, ‘31’, ‘32‘
= それぞれ非対角

テンソル成分. ただし, 対称テンソルであるため, ‘12’=’21’ を 計算していることに注意.

r(:) :real, intent(in)
: radial 方向の空間座標 [m]
z(:) :real, intent(in)
: vertical 方向の空間座標 [m]
u(size(r),size(z)) :real, intent(in)
: radial に対応する方向の 3 次元ベクトル成分
v(size(r),size(z)) :real, intent(in)
: tangential に対応する方向の 3 次元ベクトル成分
w(size(r),size(z)) :real, intent(in)
: vertical に対応する方向の 3 次元ベクトル成分
rho(size(z)) :real, intent(in)
: 水平面に平均した基本場の密度 [kg/m^3]
nuh(size(r),size(z)) :real, intent(in)
: 水平渦粘性係数
nuv(size(r),size(z)) :real, intent(in)
: 鉛直渦粘性係数
val(size(r),size(z)) :real, intent(inout)
: 計算されたテンソル成分 現在, 以下のオプションは使用していない.
undef :real, intent(in), optional
sfctau(size(r)) :real, intent(in), optional
: 地表面からのフラックス これが与えられれば, 最下層の応力はこれで置き換える.

円筒座標系におけるレイノルズ応力テンソルを計算する.

[Source]

subroutine tangent_mean_Reynolds( signal, r, z, u, v, w, rho, nuh, nuv, val, undef, sfctau )
! 円筒座標系におけるレイノルズ応力テンソルを計算する.
  use analy
  implicit none
  character(2) :: signal  ! 計算するテンソル成分.
                  ! ['11', '22', '33'] = それぞれ対角テンソル成分
                  ! ['12', '13', '21', '23', '31', '32'] = それぞれ非対角
                  ! テンソル成分. ただし, 対称テンソルであるため, '12'='21' を
                  ! 計算していることに注意.
  real, intent(in) :: r(:)  ! radial 方向の空間座標 [m]
  real, intent(in) :: z(:)  ! vertical 方向の空間座標 [m]
  real, intent(in) :: u(size(r),size(z))  ! radial に対応する方向の 3 次元ベクトル成分
  real, intent(in) :: v(size(r),size(z))  ! tangential に対応する方向の 3 次元ベクトル成分
  real, intent(in) :: w(size(r),size(z))  ! vertical に対応する方向の 3 次元ベクトル成分
  real, intent(in) :: rho(size(z))  ! 水平面に平均した基本場の密度 [kg/m^3]
  real, intent(in) :: nuh(size(r),size(z))  ! 水平渦粘性係数
  real, intent(in) :: nuv(size(r),size(z))  ! 鉛直渦粘性係数
  real, intent(inout) :: val(size(r),size(z))  ! 計算されたテンソル成分
! 現在, 以下のオプションは使用していない.
  real, intent(in), optional :: undef
  real, intent(in), optional :: sfctau(size(r))  ! 地表面からのフラックス
                 ! これが与えられれば, 最下層の応力はこれで置き換える.
  integer :: i   ! イタレーション用添字
  integer :: j   ! イタレーション用添字
  integer :: k   ! イタレーション用添字
  integer :: nr  ! 空間配列要素数 1 次元目
  integer :: nz  ! 空間配列要素数 3 次元目
  real :: dr  ! 1 次元目を微分する格子間隔 [m]
  real :: dz  ! 3 次元目を微分する格子間隔 [m]
  real, dimension(size(r),size(z)) :: tau1, tau2, tau3  ! signal 方向に
              ! 作用する 1,2,3 面に垂直な応力
  integer :: signaltau
  real :: sxx(size(r),size(z)), nu(size(r),size(z))
  real :: stau(size(r))

  dr=r(2)-r(1)
  dz=z(2)-z(1)
  nr=size(r)
  nz=size(z)

  val=0.0
  stau=0.0

  if(present(sfctau))then
     if(signal(2:2)=='3'.and.signal(1:1)/='3')then
        stau(:)=sfctau(:)
     end if
  end if

!-- [NOTE]
!-- 以下, 文字で case の or ができないため, 
!-- if 文の入れ子ではなく, if 文の並列表記で case と同じように見せかける.
!-- これはもちろん, 上から順に if をたどるが, どの場合も 2 種類以上の if に
!-- 合致しないことが既知であるために可能となる書き方であり,
!-- 並列表記した if の 2 パターン以上に合致してしまうような条件文では,
!-- case の代用には用いることができないことに注意.
!-- 本ライブラリでこのような紛らわしい表記をしている場合は必ず NOTE が入る.

  if(signal(1:2)=='12'.or.signal(1:2)=='21')then
     call tangent_mean_deform( signal, r, z, u, v, w, sxx )

     do k=1,nz
        do i=1,nr
           nu(i,k)=nuh(i,k)
        end do
     end do
  end if

  if(signal(1:2)=='23'.or.signal(1:2)=='32')then
     call tangent_mean_deform( signal, r, z, u, v, w, sxx )

     if(signal(2:2)=='3')then
        do k=1,nz
           do i=1,nr
              nu(i,k)=nuv(i,k)
           end do
        end do
     else
        do k=1,nz
           do i=1,nr
              nu(i,k)=nuh(i,k)
           end do
        end do
     end if
  end if

  if(signal(1:2)=='13'.or.signal(1:2)=='31')then
     call tangent_mean_deform( signal, r, z, u, v, w, sxx )

     if(signal(2:2)=='3')then
        do k=1,nz
           do i=1,nr
              nu(i,k)=nuv(i,k)
           end do
        end do
     else
        do k=1,nz
           do i=1,nr
              nu(i,k)=nuh(i,k)
           end do
        end do
     end if
  end if

  if(signal(1:2)=='11')then
     call tangent_mean_deform( signal, r, z, u, v, w, sxx )
     call div( r, z, u, w, val )
     do k=1,nz
        do i=1,nr
           if(r(i)/=0.0)then
              val(i,k)=val(i,k)+u(i,k)/r(i)
           end if
        end do
     end do

     do k=1,nz
        do i=1,nr
           nu(i,k)=nuh(i,k)
        end do
     end do
  end if

  if(signal(1:2)=='22')then
     call tangent_mean_deform( signal, r, z, u, v, w, sxx )
     call div( r, z, u, w, val )
     do k=1,nz
        do i=1,nr
           if(r(i)/=0.0)then
              val(i,k)=val(i,k)+u(i,k)/r(i)
           end if
        end do
     end do

     do k=1,nz
        do i=1,nr
           nu(i,k)=nuh(i,k)
        end do
     end do

  end if

  if(signal(1:2)=='33')then
     call tangent_mean_deform( signal, r, z, u, v, w, sxx )
     call div( r, z, u, w, val )
     do k=1,nz
        do i=1,nr
           if(r(i)/=0.0)then
              val(i,k)=val(i,k)+u(i,k)/r(i)
           end if
        end do
     end do

     do k=1,nz
        do i=1,nr
           nu(i,k)=nuh(i,k)
        end do
     end do

  end if

!-- 以下の式は, 最初 val = 0 で if 文の中で計算されているものとされていない
!-- ものに分かれるので, 統一の式で評価ができる.
!-- 計算されていないものについてはそもそもゼロである.

!-- 以下, 最下層は地表面フラックスを代入するかどうかのオプションのため, 別ループ

  if(present(sfctau))then
     do i=1,nr
        val(i,1)=stau(i)
     end do
  else
     do i=1,nr
        val(i,1)=rho(1)*nu(i,1)*(sxx(i,1)-(2.0/3.0)*val(i,1))
     end do
  end if

  do k=2,nz
     do i=1,nr
        val(i,k)=rho(k)*nu(i,k)*(sxx(i,k)-(2.0/3.0)*val(i,k))
     end do
  end do

end subroutine
Subroutine :
x(:) :real, intent(in)
: デカルト座標系での x 座標
y(:) :real, intent(in)
: デカルト座標系での y 座標
xc :real, intent(in)
: 接線平均する際の中心 x 成分.
yc :real, intent(in)
: 接線平均する際の中心 y 成分.
u(size(x),size(y)) :real, intent(in)
: デカルト座標系での平均化する値
r(:) :real, intent(in)
: 平均化したときの動径方向の座標(xc からの値を入れる).
theta(:) :real, intent(in)
: 平均化するときの接線方向の座標 [rad].
v(size(r),size(theta)) :real, intent(inout)
: 平均化した u のアノマリー.

任意の物理量を台風の中心から接線方向へ平均し, そのアノマリーを計算するルーチン 引数の制限として, |x(1)-xc|, |xc-x(nx)|, |y(1)-yc|, |yc-y(ny)| > r. 平均化の手順は以下のとおり. (1) nr, nt のすべての点についてそれに対応する x, y 座標値を rt_2_xy で計算. (2) その地点を含む x,y グリッドの微小領域を interpo_search_2d で検索. (3) その地点を含む 4 点が出たら, その地点でのスカラー値を 4 隅のスカラー値

    から, 重線形内挿 interpolation_2d で計算.

(4) nr x nt 個の内挿スカラー値が求まったら, nt 方向に偏差計算 Anomaly_1d 使用. 以上で各 nr について偏差値が得られる. 本ルーチンは平面極座標グリッド値での偏差計算を行う. 接線平均を行ったのち, デカルト座標の偏差へ落とし込むには, tangent_mean_anom_scal_r2c を使用.

[Source]

subroutine tangent_mean_anom_scal( x, y, xc, yc, u, r, theta, v )
  ! 任意の物理量を台風の中心から接線方向へ平均し, そのアノマリーを計算するルーチン
  ! 引数の制限として, |x(1)-xc|, |xc-x(nx)|, |y(1)-yc|, |yc-y(ny)| > r. 
  ! 平均化の手順は以下のとおり.
  ! (1) nr, nt のすべての点についてそれに対応する x, y 座標値を rt_2_xy で計算.
  ! (2) その地点を含む x,y グリッドの微小領域を interpo_search_2d で検索.
  ! (3) その地点を含む 4 点が出たら, その地点でのスカラー値を 4 隅のスカラー値
  !     から, 重線形内挿 interpolation_2d で計算.
  ! (4) nr x nt 個の内挿スカラー値が求まったら, nt 方向に偏差計算 Anomaly_1d 使用.
  ! 以上で各 nr について偏差値が得られる.
  ! 本ルーチンは平面極座標グリッド値での偏差計算を行う.
  ! 接線平均を行ったのち, デカルト座標の偏差へ落とし込むには, 
  ! tangent_mean_anom_scal_r2c を使用.
  use analy
  use max_min
  use statistics
  use Geometry
  implicit none
  real, intent(in) :: x(:)  ! デカルト座標系での x 座標
  real, intent(in) :: y(:)  ! デカルト座標系での y 座標
  real, intent(in) :: u(size(x),size(y))  ! デカルト座標系での平均化する値
  real, intent(in) :: xc  ! 接線平均する際の中心 x 成分.
  real, intent(in) :: yc  ! 接線平均する際の中心 y 成分.
  real, intent(in) :: r(:)  ! 平均化したときの動径方向の座標(xc からの値を入れる).
  real, intent(in) :: theta(:)  ! 平均化するときの接線方向の座標 [rad].
  real, intent(inout) :: v(size(r),size(theta))  ! 平均化した u のアノマリー.
  integer :: i, j, k, nx, ny, nr, nt
  real :: work(size(r),size(theta))
  real :: point(size(r),size(theta),2)
  integer :: ip(size(r),size(theta),2)
  real :: tmpx(2), tmpy(2), tmpz(2,2), inter(2)

  nx=size(x)
  ny=size(y)
  nr=size(r)
  nt=size(theta)

!-- 先の引数条件をクリアしているか確認 ---
  if(abs(x(1)-xc) < r(nr))then
     write(*,*) "error : |x(1)-xc| >= rmax. "
     stop
  else
     if(abs(x(nx)-xc) < r(nr))then
        write(*,*) "error : |x(nx)-xc| >= rmax. "
        stop
     else
        if(abs(y(1)-yc) < r(nr))then
           write(*,*) "error : |y(1)-yc| >= rmax. "
           stop
        else
           if(abs(y(ny)-yc) < r(nr))then
              write(*,*) "error : |y(ny)-yc| >= rmax. "
              stop
           end if
        end if
     end if
  end if

!-- 過程(1) ---
  do j=1,nt
     do i=1,nr
        call rt_2_xy( r(i), theta(j), point(i,j,1), point(i,j,2) )
        point(i,j,1)=xc+point(i,j,1)
        point(i,j,2)=yc+point(i,j,2)
     end do
  end do

!-- 過程(2) ---
  do j=1,nt
     do i=1,nr
        call interpo_search_2d( x, y, point(i,j,1), point(i,j,2), ip(i,j,1), ip(i,j,2) )
     end do
  end do

!-- 過程(3) ---
  do j=1,nt
     do i=1,nr
        tmpx(1)=x(ip(i,j,1))
        tmpx(2)=x(ip(i,j,1)+1)
        tmpy(1)=y(ip(i,j,2))
        tmpy(2)=y(ip(i,j,2)+1)
        tmpz(1,1)=u(ip(i,j,1),ip(i,j,2))
        tmpz(2,1)=u(ip(i,j,1)+1,ip(i,j,2))
        tmpz(1,2)=u(ip(i,j,1),ip(i,j,2)+1)
        tmpz(2,2)=u(ip(i,j,1)+1,ip(i,j,2)+1)
        inter(1)=point(i,j,1)
        inter(2)=point(i,j,2)
        call interpolation_2d( tmpx, tmpy, tmpz, inter, work(i,j) )
     end do
  end do

!-- 過程(4) ---
  do i=2,nr
     call Anomaly_1d( work(i,:), v(i,:) )
  end do
  v(1,:)=work(1,1)

end subroutine tangent_mean_anom_scal
Subroutine :
x(:) :real, intent(in)
: デカルト座標系での x 座標
y(:) :real, intent(in)
: デカルト座標系での y 座標
xc :real, intent(in)
: 接線平均する際の中心 x 成分.
yc :real, intent(in)
: 接線平均する際の中心 y 成分.
scal(size(x),size(y)) :real, intent(in)
: デカルト座標系での平均化する値.
r(:) :real, intent(in)
: 平均化したときの動径方向の座標(xc からの値を入れる).
theta(:) :real, intent(in)
: 平均化するときの接線方向の座標 [rad].
scal_anom(size(x),size(y)) :real, intent(inout)
: デカルト系でのアノマリ.
undef :real, optional
: 内挿値が見つからないときの未定義値. デフォルトでは dcl の未定義値

台風中心から接線アノマリを計算し, デカルト座標系に戻す. 接線平均ルーチンを用いて物理量を平均し, その 1 次元データを確保しておく. 同時に, デカルト系での物理量の円筒座標系での radial 位置を求める. この radial の位置における接線平均値を先の 1 次元データから内挿で求める. この求めた内挿値を元のデカルトデータから引くことでアノマリとする.

[Source]

subroutine tangent_mean_anom_scal_Cart( x, y, xc, yc, scal, r, theta, scal_anom, undef )
  ! 台風中心から接線アノマリを計算し, デカルト座標系に戻す.
  ! 接線平均ルーチンを用いて物理量を平均し, その 1 次元データを確保しておく.
  ! 同時に, デカルト系での物理量の円筒座標系での radial 位置を求める.
  ! この radial の位置における接線平均値を先の 1 次元データから内挿で求める.
  ! この求めた内挿値を元のデカルトデータから引くことでアノマリとする.
  use statistics
  implicit none
  real, intent(in) :: x(:)  ! デカルト座標系での x 座標
  real, intent(in) :: y(:)  ! デカルト座標系での y 座標
  real, intent(in) :: scal(size(x),size(y))  ! デカルト座標系での平均化する値.
  real, intent(in) :: xc  ! 接線平均する際の中心 x 成分.
  real, intent(in) :: yc  ! 接線平均する際の中心 y 成分.
  real, intent(in) :: r(:)  ! 平均化したときの動径方向の座標(xc からの値を入れる).
  real, intent(in) :: theta(:)  ! 平均化するときの接線方向の座標 [rad].
  real, intent(inout) :: scal_anom(size(x),size(y))  ! デカルト系でのアノマリ.
  real, optional :: undef  ! 内挿値が見つからないときの未定義値. デフォルトでは dcl の未定義値
  integer :: i, j, k, nx, ny, nr, nt, itmpr
  real :: tmp(size(r))
  real :: tmpr, tmp_anom, undeff

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

  nx=size(x)
  ny=size(y)
  nr=size(r)
  nt=size(theta)

  call tangent_mean_scal( x, y, xc, yc, scal, r, theta, tmp )

!-- 接線平均値を内挿し, その内挿値を引いてアノマリを求める.
  do k=1,ny
     do j=1,nx
        tmpr=sqrt((x(j)-xc)**2+(y(k)-yc)**2)
        call interpo_search_1d( r, tmpr, itmpr )
        if(nr>itmpr)then
           call interpolation_1d( r(itmpr), r(itmpr+1), tmp(itmpr), tmp(itmpr+1), tmpr, tmp_anom )
           scal_anom(j,k)=scal(j,k)-tmp_anom
        else
           scal_anom(j,k)=undeff
        end if
     end do
  end do

end subroutine tangent_mean_anom_scal_Cart
Subroutine :
charc :character(6), intent(in)
: 動径成分か接線成分かの判別, vector = 接線, scalar = 動径成分.
x(:) :real, intent(in)
: デカルト座標系での x 座標
y(:) :real, intent(in)
: デカルト座標系での y 座標
xc :real, intent(in)
: 接線平均する際の中心 x 成分.
yc :real, intent(in)
: 接線平均する際の中心 y 成分.
u1(size(x),size(y)) :real, intent(in)
: デカルト座標系での平均化する値 1
u2(size(x),size(y)) :real, intent(in)
: デカルト座標系での平均化する値 2
r(:) :real, intent(in)
: 平均化したときの動径方向の座標(xc からの値を入れる).
theta(:) :real, intent(in)
: 平均化するときの接線方向の座標 [rad].
v(size(r),size(theta)) :real, intent(inout)
: アノマリの u の値.

任意の物理量を台風の中心から接線方向へ平均アノマリを計算するルーチン 接線風速平均用. 引数の制限として, |x(1)-xc|, |xc-x(nx)|, |y(1)-yc|, |yc-y(ny)| > r. 平均化の手順は以下のとおり. (1) nr, nt のすべての点についてそれに対応する x, y 座標値を rt_2_xy で計算. (2) その地点を含む x,y グリッドの微小領域を interpo_search_2d で検索. (3) その地点を含む 4 点が出たら, その地点でのスカラー値を 4 隅のスカラー値

    から, 重線形内挿 interpolation_2d で計算.
    これと同じループ内で, 2 つのベクトル成分の平均値が得られるので,
    これらを用いて vec_prod によって中心からの位置ベクトルとの外積を計算
    . 同じループ内で中心からの距離で割って v の接線成分のみ抽出.

(4) nr x nt 個の内挿接線成分値が求まったら, nt 方向に平均計算 mean_1d 使用. 以上で各 nr について平均値が得られる.

[Source]

subroutine tangent_mean_anom_vec( charc, x, y, xc, yc, u1, u2, r, theta, v )
  ! 任意の物理量を台風の中心から接線方向へ平均アノマリを計算するルーチン
  ! 接線風速平均用.
  ! 引数の制限として, |x(1)-xc|, |xc-x(nx)|, |y(1)-yc|, |yc-y(ny)| > r. 
  ! 平均化の手順は以下のとおり.
  ! (1) nr, nt のすべての点についてそれに対応する x, y 座標値を rt_2_xy で計算.
  ! (2) その地点を含む x,y グリッドの微小領域を interpo_search_2d で検索.
  ! (3) その地点を含む 4 点が出たら, その地点でのスカラー値を 4 隅のスカラー値
  !     から, 重線形内挿 interpolation_2d で計算.
  !     これと同じループ内で, 2 つのベクトル成分の平均値が得られるので,
  !     これらを用いて vec_prod によって中心からの位置ベクトルとの外積を計算
  !     . 同じループ内で中心からの距離で割って v の接線成分のみ抽出.
  ! (4) nr x nt 個の内挿接線成分値が求まったら, nt 方向に平均計算 mean_1d 使用.
  ! 以上で各 nr について平均値が得られる.
  use analy
  use max_min
  use statistics
  use Geometry
  implicit none
  character(6), intent(in) :: charc  ! 動径成分か接線成分かの判別, vector = 接線, scalar = 動径成分.
  real, intent(in) :: x(:)  ! デカルト座標系での x 座標
  real, intent(in) :: y(:)  ! デカルト座標系での y 座標
  real, intent(in) :: u1(size(x),size(y))  ! デカルト座標系での平均化する値 1
  real, intent(in) :: u2(size(x),size(y))  ! デカルト座標系での平均化する値 2
  real, intent(in) :: xc  ! 接線平均する際の中心 x 成分.
  real, intent(in) :: yc  ! 接線平均する際の中心 y 成分.
  real, intent(in) :: r(:)  ! 平均化したときの動径方向の座標(xc からの値を入れる).
  real, intent(in) :: theta(:)  ! 平均化するときの接線方向の座標 [rad].
  real, intent(inout) :: v(size(r),size(theta))  ! アノマリの u の値.
  integer :: i, j, k, nx, ny, nr, nt
  real :: work1(size(r),size(theta),1), work2(size(r),size(theta),1), work3(size(r),size(theta),1)
  real :: posx(size(r),size(theta),1), posy(size(r),size(theta),1), posz(size(r),size(theta),1)
  real :: vecx(size(r),size(theta),1), vecy(size(r),size(theta),1), vecz(size(r),size(theta),1)
  real :: abpos(size(r),size(theta))
  real :: point(size(r),size(theta),2)
  integer :: ip(size(r),size(theta),2)
  real :: tmpx(2), tmpy(2), tmpz(2,2), inter(2)

  nx=size(x)
  ny=size(y)
  nr=size(r)
  nt=size(theta)

!-- 先の引数条件をクリアしているか確認 ---
  if(abs(x(1)-xc) < r(nr))then
     write(*,*) "error : |x(1)-xc| >= rmax. "
     stop
  else
     if(abs(x(nx)-xc) < r(nr))then
        write(*,*) "error : |x(nx)-xc| >= rmax. "
        stop
     else
        if(abs(y(1)-yc) < r(nr))then
           write(*,*) "error : |y(1)-yc| >= rmax. "
           stop
        else
           if(abs(y(ny)-yc) < r(nr))then
              write(*,*) "error : |y(ny)-yc| >= rmax. "
              stop
           end if
        end if
     end if
  end if

!-- 過程(1) ---
  do j=1,nt
     do i=1,nr
        call rt_2_xy( r(i), theta(j), point(i,j,1), point(i,j,2) )
        point(i,j,1)=xc+point(i,j,1)
        point(i,j,2)=yc+point(i,j,2)
     end do
  end do

!-- 過程(2) ---
  do j=1,nt
     do i=1,nr
        call interpo_search_2d( x, y, point(i,j,1), point(i,j,2), ip(i,j,1), ip(i,j,2) )
     end do
  end do


!-- 過程(3) ---
!-- 1. ベクトルの 2 成分について内挿した値を配列に格納 ---
  do j=1,nt
     do i=1,nr
        tmpx(1)=x(ip(i,j,1))
        tmpx(2)=x(ip(i,j,1)+1)
        tmpy(1)=y(ip(i,j,2))
        tmpy(2)=y(ip(i,j,2)+1)
        tmpz(1,1)=u1(ip(i,j,1),ip(i,j,2))
        tmpz(2,1)=u1(ip(i,j,1)+1,ip(i,j,2))
        tmpz(1,2)=u1(ip(i,j,1),ip(i,j,2)+1)
        tmpz(2,2)=u1(ip(i,j,1)+1,ip(i,j,2)+1)
        inter(1)=point(i,j,1)
        inter(2)=point(i,j,2)
        call interpolation_2d( tmpx, tmpy, tmpz, inter, work1(i,j,1) )
        tmpz(1,1)=u2(ip(i,j,1),ip(i,j,2))
        tmpz(2,1)=u2(ip(i,j,1)+1,ip(i,j,2))
        tmpz(1,2)=u2(ip(i,j,1),ip(i,j,2)+1)
        tmpz(2,2)=u2(ip(i,j,1)+1,ip(i,j,2)+1)
        call interpolation_2d( tmpx, tmpy, tmpz, inter, work2(i,j,1) )
     end do
  end do

!-- 2. 内挿した配列を使って, 内挿点での x,y 座標(その要素番号は, ix,iy に格納)
!--    の位置ベクトルとの外積を計算
  do j=1,nt
     do i=1,nr
        work3(i,j,1)=0.0
        posx(i,j,1)=point(i,j,1)-xc
        posy(i,j,1)=point(i,j,2)-yc
        posz(i,j,1)=0.0
     end do
  end do

  select case (charc)
  case ('vector')
     call vec_prod( posx, posy, posz, work1, work2, work3, vecx, vecy, vecz )
  case ('scalar')
     call dot_prod( posx, posy, posz, work1, work2, work3, vecz )
  case default
     write(*,*) "error : bad character. select 'vector', or 'scalar'."
     stop
  end select

!-- 3. ベクトルの各成分のうち, z 成分について (2 次元水平面ベクトル同士の外積)
!--    位置ベクトルの絶対値で割る. -> 接線風速成分の内挿した値が得られる.
  call abst( posx, posy, posz, abpos )

  do j=1,nt
     do i=2,nr   
        vecz(i,j,1)=vecz(i,j,1)/abpos(i,j)
     end do
  end do


!-- 過程(4) ---
  do i=2,nr
     call Anomaly_1d( vecz(i,:,1), v(i,:) )
  end do
  v(1,:)=0.0

end subroutine tangent_mean_anom_vec
Subroutine :
signal :character(2)
: 計算するテンソル成分.
‘11’, ‘22’, ‘33‘
= それぞれ対角テンソル成分
‘12’, ‘13’, ‘21’, ‘23’, ‘31’, ‘32‘
= それぞれ非対角

テンソル成分. ただし, 対称テンソルであるため, ‘12’=’21’ を 計算していることに注意.

r(:) :real, intent(in)
: radial 方向の空間座標 [m]
z(:) :real, intent(in)
: vertical 方向の空間座標 [m]
u(size(r),size(z)) :real, intent(in)
: radial に対応する方向の 3 次元ベクトル成分
v(size(r),size(z)) :real, intent(in)
: tangential に対応する方向の 3 次元ベクトル成分
w(size(r),size(z)) :real, intent(in)
: vertical に対応する方向の 3 次元ベクトル成分
val(size(r),size(z)) :real, intent(inout)
: 計算されたテンソル成分 現在, 以下のオプションは使用していない.
undef :real, intent(in), optional

デカルト座標系における変形速度テンソルを計算する.

[Source]

subroutine tangent_mean_deform( signal, r, z, u, v, w, val, undef )
! デカルト座標系における変形速度テンソルを計算する.
  use analy
  implicit none
  character(2) :: signal  ! 計算するテンソル成分.
                  ! ['11', '22', '33'] = それぞれ対角テンソル成分
                  ! ['12', '13', '21', '23', '31', '32'] = それぞれ非対角
                  ! テンソル成分. ただし, 対称テンソルであるため, '12'='21' を
                  ! 計算していることに注意.
  real, intent(in) :: r(:)  ! radial 方向の空間座標 [m]
  real, intent(in) :: z(:)  ! vertical 方向の空間座標 [m]
  real, intent(in) :: u(size(r),size(z))  ! radial に対応する方向の 3 次元ベクトル成分
  real, intent(in) :: v(size(r),size(z))  ! tangential に対応する方向の 3 次元ベクトル成分
  real, intent(in) :: w(size(r),size(z))  ! vertical に対応する方向の 3 次元ベクトル成分
  real, intent(inout) :: val(size(r),size(z))  ! 計算されたテンソル成分
! 現在, 以下のオプションは使用していない.
  real, intent(in), optional :: undef
  integer :: i   ! イタレーション用添字
  integer :: j   ! イタレーション用添字
  integer :: k   ! イタレーション用添字
  integer :: nr  ! 空間配列要素数 1 次元目
  integer :: nz  ! 空間配列要素数 2 次元目
  real :: dr  ! 1 次元目を微分する格子間隔 [m]
  real :: dz  ! 2 次元目を微分する格子間隔 [m]
  real, dimension(size(r),size(z)) :: tau1, tau2, tau3  ! signal 方向に
              ! 作用する 1,2,3 面に垂直な応力

  dr=r(2)-r(1)
  dz=z(2)-z(1)
  nr=size(r)
  nz=size(z)

!-- [NOTE]
!-- 以下, 文字で case の or ができないため, 
!-- if 文の入れ子ではなく, if 文の並列表記で case と同じように見せかける.
!-- これはもちろん, 上から順に if をたどるが, どの場合も 2 種類以上の if に
!-- 合致しないことが既知であるために可能となる書き方であり,
!-- 並列表記した if の 2 パターン以上に合致してしまうような条件文では,
!-- case の代用には用いることができないことに注意.
!-- 本ライブラリでこのような紛らわしい表記をしている場合は必ず NOTE が入る.

  if(signal(1:2)=='12'.or.signal(1:2)=='21')then
     do k=1,nz
        call grad_1d( r, v(:,k), val(:,k) )
        do i=1,nr
           if(r(i)/=0.0)then
              val(i,k)=val(i,k)-v(i,k)/r(i)
           end if
        end do
     end do
  end if

  if(signal(1:2)=='23'.or.signal(1:2)=='32')then
!$omp parallel default(shared)
!$omp do private(k)
     do k=1,nr
        call grad_1d( z, v(k,:), val(k,:) )
     end do
!$omp end do
!$omp end parallel
  end if

  if(signal(1:2)=='13'.or.signal(1:2)=='31')then
     call div( r, z, w, u, val )
  end if

  if(signal(1:2)=='11')then
!$omp parallel default(shared)
!$omp do private(k)
     do k=1,nz
        call grad_1d( r, u(:,k), val(:,k) )
        val(:,k)=2.0*val(:,k)
     end do
!$omp end do
!$omp end parallel
  end if

  if(signal(1:2)=='22')then
!$omp parallel default(shared)
!$omp do private(j,k)
     do k=1,nz
        do j=1,nr
           if(r(j)/=0.0)then
              val(j,k)=2.0*u(j,k)/r(j)
           else
              val(j,k)=0.0
           end if
        end do
     end do
!$omp end do
!$omp end parallel
  end if

  if(signal(1:2)=='33')then
!$omp parallel default(shared)
!$omp do private(j)
     do j=1,nr
        call grad_1d( z, w(j,:), val(j,:) )
        val(j,:)=2.0*val(j,:)
     end do
!$omp end do
!$omp end parallel
  end if

end subroutine
Subroutine :
x(:) :real, intent(in)
: デカルト座標系での x 座標
y(:) :real, intent(in)
: デカルト座標系での y 座標
xc :real, intent(in)
: 接線平均する際の中心 x 成分.
yc :real, intent(in)
: 接線平均する際の中心 y 成分.
u(size(x),size(y)) :real, intent(in)
: デカルト座標系での平均化する値
r(:) :real, intent(in)
: 平均化したときの動径方向の座標(xc からの値を入れる).
theta(:) :real, intent(in)
: 平均化するときの接線方向の座標 [rad].
v(size(r)) :real, intent(inout)
: 平均化した u の値.

任意の物理量を台風の中心から接線方向へ平均するルーチン このルーチンは接線風速を平均する時には用いることはできない. 接線の平均を行う際には, 別のルーチン, tangent_mean_vec の使用が必要. 引数の制限として, |x(1)-xc|, |xc-x(nx)|, |y(1)-yc|, |yc-y(ny)| > r. 平均化の手順は以下のとおり. (1) nr, nt のすべての点についてそれに対応する x, y 座標値を rt_2_xy で計算. (2) その地点を含む x,y グリッドの微小領域を interpo_search_2d で検索. (3) その地点を含む 4 点が出たら, その地点でのスカラー値を 4 隅のスカラー値

    から, 重線形内挿 interpolation_2d で計算.

(4) nr x nt 個の内挿スカラー値が求まったら, nt 方向に平均計算 mean_1d 使用. 以上で各 nr について平均値が得られる.

[Source]

subroutine tangent_mean_scal( x, y, xc, yc, u, r, theta, v )
  ! 任意の物理量を台風の中心から接線方向へ平均するルーチン
  ! このルーチンは接線風速を平均する時には用いることはできない.
  ! 接線の平均を行う際には, 別のルーチン, tangent_mean_vec の使用が必要.
  ! 引数の制限として, |x(1)-xc|, |xc-x(nx)|, |y(1)-yc|, |yc-y(ny)| > r. 
  ! 平均化の手順は以下のとおり.
  ! (1) nr, nt のすべての点についてそれに対応する x, y 座標値を rt_2_xy で計算.
  ! (2) その地点を含む x,y グリッドの微小領域を interpo_search_2d で検索.
  ! (3) その地点を含む 4 点が出たら, その地点でのスカラー値を 4 隅のスカラー値
  !     から, 重線形内挿 interpolation_2d で計算.
  ! (4) nr x nt 個の内挿スカラー値が求まったら, nt 方向に平均計算 mean_1d 使用.
  ! 以上で各 nr について平均値が得られる.
  use analy
  use max_min
  use statistics
  use Geometry
  implicit none
  real, intent(in) :: x(:)  ! デカルト座標系での x 座標
  real, intent(in) :: y(:)  ! デカルト座標系での y 座標
  real, intent(in) :: u(size(x),size(y))  ! デカルト座標系での平均化する値
  real, intent(in) :: xc  ! 接線平均する際の中心 x 成分.
  real, intent(in) :: yc  ! 接線平均する際の中心 y 成分.
  real, intent(in) :: r(:)  ! 平均化したときの動径方向の座標(xc からの値を入れる).
  real, intent(in) :: theta(:)  ! 平均化するときの接線方向の座標 [rad].
  real, intent(inout) :: v(size(r))  ! 平均化した u の値.
  integer :: i, j, k, nx, ny, nr, nt
  real :: work(size(r),size(theta))
  real :: point(size(r),size(theta),2)
  integer :: ip(size(r),size(theta),2)
  real :: tmpx(2), tmpy(2), tmpz(2,2), inter(2)

  nx=size(x)
  ny=size(y)
  nr=size(r)
  nt=size(theta)

!-- 先の引数条件をクリアしているか確認 ---
  if(abs(x(1)-xc) < r(nr))then
     write(*,*) "error : |x(1)-xc| >= rmax. "
     stop
  else
     if(abs(x(nx)-xc) < r(nr))then
        write(*,*) "error : |x(nx)-xc| >= rmax. "
        stop
     else
        if(abs(y(1)-yc) < r(nr))then
           write(*,*) "error : |y(1)-yc| >= rmax. "
           stop
        else
           if(abs(y(ny)-yc) < r(nr))then
              write(*,*) "error : |y(ny)-yc| >= rmax. "
              stop
           end if
        end if
     end if
  end if

!-- 過程(1) ---
  do j=1,nt
     do i=1,nr
        call rt_2_xy( r(i), theta(j), point(i,j,1), point(i,j,2) )
        point(i,j,1)=xc+point(i,j,1)
        point(i,j,2)=yc+point(i,j,2)
     end do
  end do

!-- 過程(2) ---
  do j=1,nt
     do i=1,nr
        call interpo_search_2d( x, y, point(i,j,1), point(i,j,2), ip(i,j,1), ip(i,j,2) )
     end do
  end do

!-- 過程(3) ---
  do j=1,nt
     do i=1,nr
        tmpx(1)=x(ip(i,j,1))
        tmpx(2)=x(ip(i,j,1)+1)
        tmpy(1)=y(ip(i,j,2))
        tmpy(2)=y(ip(i,j,2)+1)
        tmpz(1,1)=u(ip(i,j,1),ip(i,j,2))
        tmpz(2,1)=u(ip(i,j,1)+1,ip(i,j,2))
        tmpz(1,2)=u(ip(i,j,1),ip(i,j,2)+1)
        tmpz(2,2)=u(ip(i,j,1)+1,ip(i,j,2)+1)
        inter(1)=point(i,j,1)
        inter(2)=point(i,j,2)
        call interpolation_2d( tmpx, tmpy, tmpz, inter, work(i,j) )
     end do
  end do

!-- 過程(4) ---
  do i=2,nr
     call Mean_1d( work(i,:), v(i) )
  end do
  v(1)=work(1,1)

end subroutine tangent_mean_scal
Subroutine :
signal :character(1)
: 円筒座標系の何番目の乱流成分かを判定する.
1
= 円筒座標における radial 座標成分 (方程式 vr 成分)
2
= 円筒座標における tangential 座標成分 (方程式 vt 成分)
3
= 円筒座標における vertical 座標成分 (方程式 w 成分)
r(:) :real, intent(in)
: 動径方向の位置座標 [m]
z(:) :real, intent(in)
: 鉛直方向の位置座標 [m]
u(size(r),size(z)) :real, intent(in)
: x に対応する方向の 2 次元ベクトル成分
v(size(r),size(z)) :real, intent(in)
: y に対応する方向の 2 次元ベクトル成分
w(size(r),size(z)) :real, intent(in)
: y に対応する方向の 2 次元ベ>クトル成分
rho(size(z)) :real, intent(in)
: 水平面に平均した基本場の密度 [kg/m^3]
nuh(size(r),size(z)) :real, intent(in)
: 水平渦粘性係数
nuv(size(r),size(z)) :real, intent(in)
: 鉛直渦粘性係数
val(size(r),size(z)) :real, intent(inout)
: 乱流フラックス
undef :real, intent(in), optional
sfctau(size(r)) :real, intent(in), optional
: 地表面からのフラックス これが与えられれば, 最下層の応力はこれで置き換える.
 接線平均した乱流フラックスを計算する.
 接線平均しているので, tau_{*2} 成分 (\theta 微分の成分) は含まれない.

[Source]

subroutine tangent_mean_turb( signal, r, z, u, v, w, rho, nuh, nuv, val, undef, sfctau )
!  接線平均した乱流フラックスを計算する.
!  接線平均しているので, tau_{*2} 成分 (\theta 微分の成分) は含まれない.
  use analy
  implicit none
  character(1) :: signal  ! 円筒座標系の何番目の乱流成分かを判定する.
                  ! [1] = 円筒座標における radial 座標成分 (方程式 vr 成分)
                  ! [2] = 円筒座標における tangential 座標成分 (方程式 vt 成分)
                  ! [3] = 円筒座標における vertical 座標成分 (方程式 w 成分)
  real, intent(in) :: r(:)  ! 動径方向の位置座標 [m]
  real, intent(in) :: z(:)  ! 鉛直方向の位置座標 [m]
  real, intent(in) :: u(size(r),size(z))  ! x に対応する方向の 2 次元ベクトル成分
  real, intent(in) :: v(size(r),size(z))  ! y に対応する方向の 2 次元ベクトル成分
  real, intent(in) :: w(size(r),size(z))  ! y に対応する方向の 2 次元ベ>クトル成分
  real, intent(in) :: rho(size(z))  ! 水平面に平均した基本場の密度 [kg/m^3]
  real, intent(in) :: nuh(size(r),size(z))  ! 水平渦粘性係数
  real, intent(in) :: nuv(size(r),size(z))  ! 鉛直渦粘性係数
  real, intent(inout) :: val(size(r),size(z))  ! 乱流フラックス
  real, intent(in), optional :: undef
  real, intent(in), optional :: sfctau(size(r))  ! 地表面からのフラックス
                 ! これが与えられれば, 最下層の応力はこれで置き換える.
  integer :: i   ! イタレーション用添字
  integer :: j   ! イタレーション用添字
  integer :: k   ! イタレーション用添字
  integer :: id   ! イタレーション用添字
  integer :: nr  ! 空間配列要素数 1 次元目
  integer :: nz  ! 空間配列要素数 2 次元目
  real :: dr  ! 1 次元目を微分する格子間隔 [m]
  real :: dz  ! 2 次元目を微分する格子間隔 [m]
  character(1) :: signaltau(3)
  real, dimension(size(r),size(z),3) :: tau  ! signal 方向に
              ! 作用する 1,2,3 面に垂直な応力
  real, dimension(size(r),size(z)) :: tmp
  real, dimension(size(r)) :: stau

  signaltau=(/ '1', '2', '3' /)

  dr=r(2)-r(1)
  dz=z(2)-z(1)
  nr=size(r)
  nz=size(z) 

  val=0.0

  do id=1,3
     if(id/=2)then  ! tau_{*2} 成分はゼロなので, 計算しない.
        if(present(sfctau))then
           stau(:)=sfctau(:)
           call tangent_mean_Reynolds( signal//signaltau(id), r, z, u, v, w, rho, nuh, nuv, tau(:,:,id), sfctau=stau )
        else
           call tangent_mean_Reynolds( signal//signaltau(id), r, z, u, v, w, rho, nuh, nuv, tau(:,:,id) )
        end if
     end if
  end do

!-- (signal, 1) 成分の計算
  do k=1,nz
     call grad_1d( r, tau(:,k,1), tmp(:,k))
     do i=1,nr
        if(r(i)/=0.0)then
           val(i,k)=tmp(i,k)+val(i,k)+tau(i,k,1)/r(i)
        else
           val(i,k)=tmp(i,k)+val(i,k)
        end if
     end do
  end do

!-- (signal, 3) 成分の計算
  do i=1,nr
     call grad_1d( z, tau(i,:,3), tmp(i,:))
     do k=1,nz
        val(i,k)=tmp(i,k)+val(i,k)
     end do
  end do



end subroutine
Subroutine :
charc :character(6), intent(in)
: 動径成分か接線成分かの判別, vector = 接線, scalar = 動径成分.
x(:) :real, intent(in)
: デカルト座標系での x 座標
y(:) :real, intent(in)
: デカルト座標系での y 座標
xc :real, intent(in)
: 接線平均する際の中心 x 成分.
yc :real, intent(in)
: 接線平均する際の中心 y 成分.
u1(size(x),size(y)) :real, intent(in)
: デカルト座標系での平均化する値 1
u2(size(x),size(y)) :real, intent(in)
: デカルト座標系での平均化する値 2
r(:) :real, intent(in)
: 平均化したときの動径方向の座標(xc からの値を入れる).
theta(:) :real, intent(in)
: 平均化するときの接線方向の座標 [rad].
v(size(r)) :real, intent(inout)
: 平均化した u の値.

任意の物理量を台風の中心から接線方向へ平均するルーチン 接線風速平均用. 引数の制限として, |x(1)-xc|, |xc-x(nx)|, |y(1)-yc|, |yc-y(ny)| > r. 平均化の手順は以下のとおり. (1) nr, nt のすべての点についてそれに対応する x, y 座標値を rt_2_xy で計算. (2) その地点を含む x,y グリッドの微小領域を interpo_search_2d で検索. (3) その地点を含む 4 点が出たら, その地点でのスカラー値を 4 隅のスカラー値

    から, 重線形内挿 interpolation_2d で計算.
    これと同じループ内で, 2 つのベクトル成分の平均値が得られるので,
    これらを用いて vec_prod によって中心からの位置ベクトルとの外積を計算
    . 同じループ内で中心からの距離で割って v の接線成分のみ抽出.

(4) nr x nt 個の内挿接線成分値が求まったら, nt 方向に平均計算 mean_1d 使用. 以上で各 nr について平均値が得られる.

[Source]

subroutine tangent_mean_vec( charc, x, y, xc, yc, u1, u2, r, theta, v )
  ! 任意の物理量を台風の中心から接線方向へ平均するルーチン
  ! 接線風速平均用.
  ! 引数の制限として, |x(1)-xc|, |xc-x(nx)|, |y(1)-yc|, |yc-y(ny)| > r. 
  ! 平均化の手順は以下のとおり.
  ! (1) nr, nt のすべての点についてそれに対応する x, y 座標値を rt_2_xy で計算.
  ! (2) その地点を含む x,y グリッドの微小領域を interpo_search_2d で検索.
  ! (3) その地点を含む 4 点が出たら, その地点でのスカラー値を 4 隅のスカラー値
  !     から, 重線形内挿 interpolation_2d で計算.
  !     これと同じループ内で, 2 つのベクトル成分の平均値が得られるので,
  !     これらを用いて vec_prod によって中心からの位置ベクトルとの外積を計算
  !     . 同じループ内で中心からの距離で割って v の接線成分のみ抽出.
  ! (4) nr x nt 個の内挿接線成分値が求まったら, nt 方向に平均計算 mean_1d 使用.
  ! 以上で各 nr について平均値が得られる.
  use analy
  use max_min
  use statistics
  use Geometry
  implicit none
  character(6), intent(in) :: charc  ! 動径成分か接線成分かの判別, vector = 接線, scalar = 動径成分.
  real, intent(in) :: x(:)  ! デカルト座標系での x 座標
  real, intent(in) :: y(:)  ! デカルト座標系での y 座標
  real, intent(in) :: u1(size(x),size(y))  ! デカルト座標系での平均化する値 1
  real, intent(in) :: u2(size(x),size(y))  ! デカルト座標系での平均化する値 2
  real, intent(in) :: xc  ! 接線平均する際の中心 x 成分.
  real, intent(in) :: yc  ! 接線平均する際の中心 y 成分.
  real, intent(in) :: r(:)  ! 平均化したときの動径方向の座標(xc からの値を入れる).
  real, intent(in) :: theta(:)  ! 平均化するときの接線方向の座標 [rad].
  real, intent(inout) :: v(size(r))  ! 平均化した u の値.
  integer :: i, j, k, nx, ny, nr, nt
  real :: work1(size(r),size(theta),1), work2(size(r),size(theta),1), work3(size(r),size(theta),1)
  real :: posx(size(r),size(theta),1), posy(size(r),size(theta),1), posz(size(r),size(theta),1)
  real :: vecx(size(r),size(theta),1), vecy(size(r),size(theta),1), vecz(size(r),size(theta),1)
  real :: abpos(size(r),size(theta))
  real :: point(size(r),size(theta),2)
  integer :: ip(size(r),size(theta),2)
  real :: tmpx(2), tmpy(2), tmpz(2,2), inter(2)

  nx=size(x)
  ny=size(y)
  nr=size(r)
  nt=size(theta)

!-- 先の引数条件をクリアしているか確認 ---
  if(abs(x(1)-xc) < r(nr))then
     write(*,*) "error : |x(1)-xc| >= rmax. "
     stop
  else
     if(abs(x(nx)-xc) < r(nr))then
        write(*,*) "error : |x(nx)-xc| >= rmax. "
        stop
     else
        if(abs(y(1)-yc) < r(nr))then
           write(*,*) "error : |y(1)-yc| >= rmax. "
           stop
        else
           if(abs(y(ny)-yc) < r(nr))then
              write(*,*) "error : |y(ny)-yc| >= rmax. "
              stop
           end if
        end if
     end if
  end if

!-- 過程(1) ---
  do j=1,nt
     do i=1,nr
        call rt_2_xy( r(i), theta(j), point(i,j,1), point(i,j,2) )
        point(i,j,1)=xc+point(i,j,1)
        point(i,j,2)=yc+point(i,j,2)
     end do
  end do

!-- 過程(2) ---
  do j=1,nt
     do i=1,nr
        call interpo_search_2d( x, y, point(i,j,1), point(i,j,2), ip(i,j,1), ip(i,j,2) )
     end do
  end do


!-- 過程(3) ---
!-- 1. ベクトルの 2 成分について内挿した値を配列に格納 ---
  do j=1,nt
     do i=1,nr
        tmpx(1)=x(ip(i,j,1))
        tmpx(2)=x(ip(i,j,1)+1)
        tmpy(1)=y(ip(i,j,2))
        tmpy(2)=y(ip(i,j,2)+1)
        tmpz(1,1)=u1(ip(i,j,1),ip(i,j,2))
        tmpz(2,1)=u1(ip(i,j,1)+1,ip(i,j,2))
        tmpz(1,2)=u1(ip(i,j,1),ip(i,j,2)+1)
        tmpz(2,2)=u1(ip(i,j,1)+1,ip(i,j,2)+1)
        inter(1)=point(i,j,1)
        inter(2)=point(i,j,2)
        call interpolation_2d( tmpx, tmpy, tmpz, inter, work1(i,j,1) )
        tmpz(1,1)=u2(ip(i,j,1),ip(i,j,2))
        tmpz(2,1)=u2(ip(i,j,1)+1,ip(i,j,2))
        tmpz(1,2)=u2(ip(i,j,1),ip(i,j,2)+1)
        tmpz(2,2)=u2(ip(i,j,1)+1,ip(i,j,2)+1)
        call interpolation_2d( tmpx, tmpy, tmpz, inter, work2(i,j,1) )
     end do
  end do

!-- 2. 内挿した配列を使って, 内挿点での x,y 座標(その要素番号は, ix,iy に格納)
!--    の位置ベクトルとの外積を計算
  do j=1,nt
     do i=1,nr
        work3(i,j,1)=0.0
        posx(i,j,1)=point(i,j,1)-xc
        posy(i,j,1)=point(i,j,2)-yc
        posz(i,j,1)=0.0
     end do
  end do

  select case (charc)
  case ('vector')
     call vec_prod( posx, posy, posz, work1, work2, work3, vecx, vecy, vecz )
  case ('scalar')
     call dot_prod( posx, posy, posz, work1, work2, work3, vecz )
  case default
     write(*,*) "error : bad character. select 'vector', or 'scalar'."
     stop
  end select

!-- 3. ベクトルの各成分のうち, z 成分について (2 次元水平面ベクトル同士の外積)
!--    位置ベクトルの絶対値で割る. -> 接線風速成分の内挿した値が得られる.
  call abst( posx, posy, posz, abpos )

  do j=1,nt
     do i=2,nr   
        vecz(i,j,1)=vecz(i,j,1)/abpos(i,j)
     end do
  end do


!-- 過程(4) ---
  do i=2,nr
     call Mean_1d( vecz(i,:,1), v(i) )
  end do
  v(1)=0.0

end subroutine tangent_mean_vec