Class Trajectory
In: trajectory.f90

トラジェクトリ解析を行うルーチン

Methods

Included Modules

Statistics

Public Instance methods

Subroutine :
dt :real, intent(in)
: 計算する時間間隔 [s]
step :integer, intent(in)
: 計算を行うステップ数
ini_x :real, intent(in)
: 初期位置 x [m]
ini_y :real, intent(in)
: 初期位置 y [m]
x(:) :real, intent(in)
: x 方向の座標値
y(:) :real, intent(in)
: y 方向の座標値
u(size(x),size(y)) :real, intent(in)
: x 方向の風速成分 [m/s]
v(size(x),size(y)) :real, intent(in)
: x 方向の風速成分 [m/s]
trajx(step) :real, intent(inout)
: トラジェクトリの位置座標 x [m]
trajy(step) :real, intent(inout)
: トラジェクトリの位置座標 y [m]
opt :character(*), intent(in), optional
: 時間積分スキーム
undef :real, intent(in), optional
: 未定義

前方流跡線解析ルーチン 時間積分スキームについては, デフォルトではオイラー法で計算. ‘HO1’= 修正オイラー法(ホインスキーム) ‘ME1’= 改良オイラー法 ‘RK3’= 3 次ルンゲクッタ ‘RK4’= 4 次ルンゲクッタ

[Source]

subroutine Forward_Traject_2d( dt, step, ini_x, ini_y, x, y, u, v, trajx, trajy, opt, undef )
! 前方流跡線解析ルーチン
! 時間積分スキームについては, デフォルトではオイラー法で計算.
! 'HO1'= 修正オイラー法(ホインスキーム)
! 'ME1'= 改良オイラー法
! 'RK3'= 3 次ルンゲクッタ
! 'RK4'= 4 次ルンゲクッタ
  implicit none
  real, intent(in) :: dt  ! 計算する時間間隔 [s]
  integer, intent(in) :: step  ! 計算を行うステップ数
  real, intent(in) :: ini_x  ! 初期位置 x [m]
  real, intent(in) :: ini_y  ! 初期位置 y [m]
  real, intent(in) :: x(:)  ! x 方向の座標値
  real, intent(in) :: y(:)  ! y 方向の座標値
  real, intent(in) :: u(size(x),size(y))  ! x 方向の風速成分 [m/s]
  real, intent(in) :: v(size(x),size(y))  ! x 方向の風速成分 [m/s]
  real, intent(inout) :: trajx(step)  ! トラジェクトリの位置座標 x [m]
  real, intent(inout) :: trajy(step)  ! トラジェクトリの位置座標 y [m]
  character(*), intent(in), optional :: opt  ! 時間積分スキーム
  real, intent(in), optional :: undef  ! 未定義
  character(3) :: option
  real :: defun
!-- 以下, 作業用配列の定義
  integer :: i, j, k, l, m, n, id
  integer :: nx, ny
  real :: k1, k2, k3, k4, l1, l2, l3, l4
  real :: u1, v1, x1, y1
  real :: inter_p(2)  ! 内挿点での位置座標(毎回更新)
  real :: inter_v(2)  ! 内挿点での速度の値(1=u, 2=v)

!-- 汎用性を持たせるために, x, y がどのような範囲をとっても計算可能な仕様にする
!-- そのため, x, y の範囲の要素を限定する
!  call (1)
!  call (nx)
!  call (1)
!  call (ny)

  nx=size(x)
  ny=size(y)

  trajx(1)=ini_x
  trajy(1)=ini_y

  if(present(opt))then
     option=opt
  else
     option='EU1'
  end if

  if(present(undef))then
     defun=undef
  else
     defun=0.0
  end if

  select case (option)
  case ('EU1')
     do i=1,step-1
!-- 初期位置における風速を内挿
!-- 初期位置の定義
        inter_p=(/trajx(i), trajy(i)/)
!-- 初期位置の近傍 4 点を計算.
        call interpo_search_2d( x, y, trajx(i), trajy(i), m, n )  
!-- u 成分についての重線形内挿
        call interpolation_2d( x(m:m+1), y(n:n+1), u(m:m+1,n:n+1), inter_p, inter_v(1) )  
!-- v 成分についての重線形内挿
        call interpolation_2d( x(m:m+1), y(n:n+1), v(m:m+1,n:n+1), inter_p, inter_v(2) )  

        x1=trajx(i)+dt*inter_v(1)
        y1=trajy(i)+dt*inter_v(2)
!-- 計算した traj が領域内に存在しているか確認
        if(x1<x(1).or.x1>x(nx).or.y1<y(1).or.y1>y(ny))then
           write(*,*) "****** ERROR ******"
           write(*,*) "This point does not exist in the region."
           write(*,*) "Stop.!!"
!-- 計算領域外に出た場合は, 以降のステップにすべて defun を代入
           do id=i,step-1
              trajx(id+1)=defun
              trajy(id+1)=defun
           end do
           exit
        end if
        trajx(i+1)=x1
        trajy(i+1)=y1
     end do

  case ('HO1')
     do i=1,step-1
        inter_p=(/trajx(i), trajy(i)/)
        call interpo_search_2d( x, y, trajx(i), trajy(i), m, n )  
        call interpolation_2d( x(m:m+1), y(n:n+1), u(m:m+1,n:n+1), inter_p, inter_v(1) )  
        call interpolation_2d( x(m:m+1), y(n:n+1), v(m:m+1,n:n+1), inter_p, inter_v(2) )  

        k1=dt*inter_v(1)
        l1=dt*inter_v(2)
!-- 一時的な流跡線の位置を計算
        x1=trajx(i)+k1
        y1=trajy(i)+l1
!-- この点での速度場を内挿
        inter_p=(/x1, y1/)
        call interpo_search_2d( x, y, x1, y1, m, n )  
        call interpolation_2d( x(m:m+1), y(n:n+1), u(m:m+1,n:n+1), inter_p, inter_v(1) )  
        call interpolation_2d( x(m:m+1), y(n:n+1), v(m:m+1,n:n+1), inter_p, inter_v(2) )  

        k2=dt*inter_v(1)
        l2=dt*inter_v(2)
        x1=trajx(i)+0.5*(k1+k2)
        y1=trajy(i)+0.5*(l1+l2)

!-- 計算した traj が領域内に存在しているか確認
        if(x1<x(1).or.x1>x(nx).or.y1<y(1).or.y1>y(ny))then
           write(*,*) "****** ERROR ******"
           write(*,*) "This point does not exist in the region."
           write(*,*) "Stop.!!"
!-- 計算領域外に出た場合は, 以降のステップにすべて defun を代入
           do id=i,step-1
              trajx(id+1)=defun
              trajy(id+1)=defun
           end do
           exit
        end if
        trajx(i+1)=x1
        trajy(i+1)=y1
     end do


  case ('ME1')
     do i=1,step-1
        inter_p=(/trajx(i), trajy(i)/)
        call interpo_search_2d( x, y, trajx(i), trajy(i), m, n )  
        call interpolation_2d( x(m:m+1), y(n:n+1), u(m:m+1,n:n+1), inter_p, inter_v(1) )  
        call interpolation_2d( x(m:m+1), y(n:n+1), v(m:m+1,n:n+1), inter_p, inter_v(2) )  

        k1=dt*inter_v(1)
        l1=dt*inter_v(2)
!-- 一時的な流跡線の位置を計算
        x1=trajx(i)+0.5*k1
        y1=trajy(i)+0.5*l1
!-- この点での速度場を内挿
        inter_p=(/x1, y1/)
        call interpo_search_2d( x, y, x1, y1, m, n )  
        call interpolation_2d( x(m:m+1), y(n:n+1), u(m:m+1,n:n+1), inter_p, inter_v(1) )  
        call interpolation_2d( x(m:m+1), y(n:n+1), v(m:m+1,n:n+1), inter_p, inter_v(2) )  

        k2=dt*inter_v(1)
        l2=dt*inter_v(2)
        x1=trajx(i)+k2
        y1=trajy(i)+l2

!-- 計算した traj が領域内に存在しているか確認
        if(x1<x(1).or.x1>x(nx).or.y1<y(1).or.y1>y(ny))then
           write(*,*) "****** ERROR ******"
           write(*,*) "This point does not exist in the region."
           write(*,*) "Stop.!!"
!-- 計算領域外に出た場合は, 以降のステップにすべて defun を代入
           do id=i,step-1
              trajx(id+1)=defun
              trajy(id+1)=defun
           end do
           exit
        end if
        trajx(i+1)=x1
        trajy(i+1)=y1
     end do



  case ('RK3')
     do i=1,step-1
        inter_p=(/trajx(i), trajy(i)/)
        call interpo_search_2d( x, y, trajx(i), trajy(i), m, n )  
        call interpolation_2d( x(m:m+1), y(n:n+1), u(m:m+1,n:n+1), inter_p, inter_v(1) )  
        call interpolation_2d( x(m:m+1), y(n:n+1), v(m:m+1,n:n+1), inter_p, inter_v(2) )  

        k1=dt*inter_v(1)
        l1=dt*inter_v(2)
!-- 一時的な流跡線の位置を計算
        x1=trajx(i)+0.5*k1
        y1=trajy(i)+0.5*l1
!-- この点での速度場を内挿
        inter_p=(/x1, y1/)
        call interpo_search_2d( x, y, x1, y1, m, n )  
        call interpolation_2d( x(m:m+1), y(n:n+1), u(m:m+1,n:n+1), inter_p, inter_v(1) )  
        call interpolation_2d( x(m:m+1), y(n:n+1), v(m:m+1,n:n+1), inter_p, inter_v(2) )  

        k2=dt*inter_v(1)
        l2=dt*inter_v(2)
        x1=trajx(i)-k1+2.0*k2
        y1=trajy(i)-l1+2.0*l2

!-- この点での速度場を内挿
        inter_p=(/x1, y1/)
        call interpo_search_2d( x, y, x1, y1, m, n )  
        call interpolation_2d( x(m:m+1), y(n:n+1), u(m:m+1,n:n+1), inter_p, inter_v(1) )  
        call interpolation_2d( x(m:m+1), y(n:n+1), v(m:m+1,n:n+1), inter_p, inter_v(2) )  

        k3=dt*inter_v(1)
        l3=dt*inter_v(2)
        x1=trajx(i)+(1.0/6.0)*(k1+4.0*k2+k3)
        y1=trajy(i)+(1.0/6.0)*(l1+4.0*l2+l3)

!-- 計算した traj が領域内に存在しているか確認
        if(x1<x(1).or.x1>x(nx).or.y1<y(1).or.y1>y(ny))then
           write(*,*) "****** ERROR ******"
           write(*,*) "This point does not exist in the region."
           write(*,*) "Stop.!!"
!-- 計算領域外に出た場合は, 以降のステップにすべて defun を代入
           do id=i,step-1
              trajx(id+1)=defun
              trajy(id+1)=defun
           end do
           exit
        end if
        trajx(i+1)=x1
        trajy(i+1)=y1
     end do


  case ('RK4')
     do i=1,step-1
        inter_p=(/trajx(i), trajy(i)/)
        call interpo_search_2d( x, y, trajx(i), trajy(i), m, n )  
        call interpolation_2d( x(m:m+1), y(n:n+1), u(m:m+1,n:n+1), inter_p, inter_v(1) )  
        call interpolation_2d( x(m:m+1), y(n:n+1), v(m:m+1,n:n+1), inter_p, inter_v(2) )  

        k1=dt*inter_v(1)
        l1=dt*inter_v(2)
!-- 一時的な流跡線の位置を計算
        x1=trajx(i)+0.5*k1
        y1=trajy(i)+0.5*l1
!-- この点での速度場を内挿
        inter_p=(/x1, y1/)
        call interpo_search_2d( x, y, x1, y1, m, n )  
        call interpolation_2d( x(m:m+1), y(n:n+1), u(m:m+1,n:n+1), inter_p, inter_v(1) )  
        call interpolation_2d( x(m:m+1), y(n:n+1), v(m:m+1,n:n+1), inter_p, inter_v(2) )  

        k2=dt*inter_v(1)
        l2=dt*inter_v(2)
        x1=trajx(i)+0.5*k2
        y1=trajy(i)+0.5*l2

!-- この点での速度場を内挿
        inter_p=(/x1, y1/)
        call interpo_search_2d( x, y, x1, y1, m, n )  
        call interpolation_2d( x(m:m+1), y(n:n+1), u(m:m+1,n:n+1), inter_p, inter_v(1) )  
        call interpolation_2d( x(m:m+1), y(n:n+1), v(m:m+1,n:n+1), inter_p, inter_v(2) )  

        k3=dt*inter_v(1)
        l3=dt*inter_v(2)

        x1=trajx(i)+k3
        y1=trajy(i)+l3

!-- この点での速度場を内挿
        inter_p=(/x1, y1/)
        call interpo_search_2d( x, y, x1, y1, m, n )  
        call interpolation_2d( x(m:m+1), y(n:n+1), u(m:m+1,n:n+1), inter_p, inter_v(1) )  
        call interpolation_2d( x(m:m+1), y(n:n+1), v(m:m+1,n:n+1), inter_p, inter_v(2) )  

        k4=dt*inter_v(1)
        l4=dt*inter_v(2)
        x1=trajx(i)+(1.0/6.0)*(k1+2.0*k2+2.0*k3+k4)
        y1=trajy(i)+(1.0/6.0)*(l1+2.0*l2+2.0*l3+l4)

!-- 計算した traj が領域内に存在しているか確認
        if(x1<x(1).or.x1>x(nx).or.y1<y(1).or.y1>y(ny))then
           write(*,*) "****** ERROR ******"
           write(*,*) "This point does not exist in the region."
           write(*,*) "Stop.!!"
!-- 計算領域外に出た場合は, 以降のステップにすべて defun を代入
           do id=i,step-1
              trajx(id+1)=defun
              trajy(id+1)=defun
           end do
           exit
        end if
        trajx(i+1)=x1
        trajy(i+1)=y1
     end do



  end select

end subroutine Forward_Traject_2d
Subroutine :
dt :real, intent(in)
: 計算する時間間隔 [s]
step :integer, intent(in)
: 計算を行うステップ数
ini_x :real, intent(in)
: 初期位置 x [m]
ini_y :real, intent(in)
: 初期位置 y [m]
ini_z :real, intent(in)
: 初期位置 z [m]
x(:) :real, intent(in)
: x 方向の座標値
y(:) :real, intent(in)
: y 方向の座標値
z(:) :real, intent(in)
: y 方向の座標値
u(size(x),size(y),size(z)) :real, intent(in)
: x 方向の風速成分 [m/s]
v(size(x),size(y),size(z)) :real, intent(in)
: x 方向の風速成分 [m/s]
w(size(x),size(y),size(z)) :real, intent(in)
: z 方向の風速成分 [m/s]
trajx(step) :real, intent(inout)
: トラジェクトリの位置座標 x [m]
trajy(step) :real, intent(inout)
: トラジェクトリの位置座標 y [m]
trajz(step) :real, intent(inout)
: トラジェクトリの位置座標 z [m]
opt :character(*), intent(in), optional
: 時間積分スキーム

前方流跡線解析ルーチン 時間積分スキームについては, デフォルトではオイラー法で計算. ‘HO1’= 修正オイラー法(ホインスキーム) ‘ME1’= 改良オイラー法 ‘RK3’= 3 次ルンゲクッタ ‘RK4’= 4 次ルンゲクッタ

[Source]

subroutine Forward_Traject_3d( dt, step, ini_x, ini_y, ini_z, x, y, z, u, v, w, trajx, trajy, trajz, opt )
! 前方流跡線解析ルーチン
! 時間積分スキームについては, デフォルトではオイラー法で計算.
! 'HO1'= 修正オイラー法(ホインスキーム)
! 'ME1'= 改良オイラー法
! 'RK3'= 3 次ルンゲクッタ
! 'RK4'= 4 次ルンゲクッタ
  implicit none
  real, intent(in) :: dt  ! 計算する時間間隔 [s]
  integer, intent(in) :: step  ! 計算を行うステップ数
  real, intent(in) :: ini_x  ! 初期位置 x [m]
  real, intent(in) :: ini_y  ! 初期位置 y [m]
  real, intent(in) :: ini_z  ! 初期位置 z [m]
  real, intent(in) :: x(:)  ! x 方向の座標値
  real, intent(in) :: y(:)  ! y 方向の座標値
  real, intent(in) :: z(:)  ! y 方向の座標値
  real, intent(in) :: u(size(x),size(y),size(z))  ! x 方向の風速成分 [m/s]
  real, intent(in) :: v(size(x),size(y),size(z))  ! x 方向の風速成分 [m/s]
  real, intent(in) :: w(size(x),size(y),size(z))  ! z 方向の風速成分 [m/s]
  real, intent(inout) :: trajx(step)  ! トラジェクトリの位置座標 x [m]
  real, intent(inout) :: trajy(step)  ! トラジェクトリの位置座標 y [m]
  real, intent(inout) :: trajz(step)  ! トラジェクトリの位置座標 z [m]
  character(*), intent(in), optional :: opt  ! 時間積分スキーム
  character(3) :: option
  real :: defun
!-- 以下, 作業用配列の定義
  integer :: i, j, k, l, m, n, nx, ny, nz, id
  real :: k1, k2, k3, k4, l1, l2, l3, l4
  real :: u1, v1, w1, x1, y1, z1, n1, n2, n3, n4
  real :: inter_p(3)  ! 内挿点での位置座標(毎回更新)
  real :: inter_v(3)  ! 内挿点での速度の値(1=u, 2=v)

!-- 汎用性を持たせるために, x, y がどのような範囲をとっても計算可能な仕様にする
!-- そのため, x, y の範囲の要素を限定する
!  call (1)
!  call (nx)
!  call (1)
!  call (ny)

  nx=size(x)
  ny=size(y)
  nz=size(z)

  trajx(1)=ini_x
  trajy(1)=ini_y
  trajz(1)=ini_z

  if(present(opt))then
     option=opt
  else
     option='EU1'
  end if

  select case (option)
  case ('EU1')
     do i=1,step-1
!-- 初期位置における風速を内挿
!-- 初期位置の定義
        inter_p=(/trajx(i), trajy(i), trajz(i)/)
!-- 初期位置の近傍 4 点を計算.
        call interpo_search_3d( x, y, z, trajx(i), trajy(i), trajz(i), m, n, l )  
!-- u 成分についての重線形内挿
        call interpolation_3d( x(m:m+1), y(n:n+1), z(l:l+1), u(m:m+1,n:n+1,l:l+1), inter_p, inter_v(1) )  
!-- v 成分についての重線形内挿
        call interpolation_3d( x(m:m+1), y(n:n+1), z(l:l+1), v(m:m+1,n:n+1,l:l+1), inter_p, inter_v(2) )  
!-- w 成分についての重線形内挿
        call interpolation_3d( x(m:m+1), y(n:n+1), z(l:l+1), w(m:m+1,n:n+1,l:l+1), inter_p, inter_v(3) )  

        x1=trajx(i)+dt*inter_v(1)
        y1=trajy(i)+dt*inter_v(2)
        z1=trajz(i)+dt*inter_v(3)
!-- 計算した traj が領域内に存在しているか確認
        if(x1<x(1).or.x1>x(nx).or.y1<y(1).or.y1>y(ny).or.z1<z(1).or.z1>z(nz))then
           write(*,*) "****** ERROR ******"
           write(*,*) "This point does not exist in the region."
           write(*,*) "Stop.!!"
!-- 計算領域外に出た場合は, 以降のステップにすべて defun を代入
           do id=i,step-1
              trajx(id+1)=defun
              trajy(id+1)=defun
              trajz(id+1)=defun
           end do
           exit
        end if
        trajx(i+1)=x1
        trajy(i+1)=y1
        trajz(i+1)=z1
     end do

  case ('HO1')
     do i=1,step-1
        inter_p=(/trajx(i), trajy(i), trajz(i)/)
        call interpo_search_3d( x, y, z, trajx(i), trajy(i), trajz(i), m, n, l )  
        call interpolation_3d( x(m:m+1), y(n:n+1), z(l:l+1), u(m:m+1,n:n+1,l:l+1), inter_p, inter_v(1) )  
        call interpolation_3d( x(m:m+1), y(n:n+1), z(l:l+1), v(m:m+1,n:n+1,l:l+1), inter_p, inter_v(2) )  
        call interpolation_3d( x(m:m+1), y(n:n+1), z(l:l+1), w(m:m+1,n:n+1,l:l+1), inter_p, inter_v(3) )  

        k1=dt*inter_v(1)
        l1=dt*inter_v(2)
        n1=dt*inter_v(3)
!-- 一時的な流跡線の位置を計算
        x1=trajx(i)+k1
        y1=trajy(i)+l1
        z1=trajz(i)+n1
!-- この点での速度場を内挿
        inter_p=(/x1, y1, z1/)
        call interpo_search_3d( x, y, z, x1, y1, z1, m, n, l )  
        call interpolation_3d( x(m:m+1), y(n:n+1), z(l:l+1), u(m:m+1,n:n+1,l:l+1), inter_p, inter_v(1) )  
        call interpolation_3d( x(m:m+1), y(n:n+1), z(l:l+1), v(m:m+1,n:n+1,l:l+1), inter_p, inter_v(2) )  
        call interpolation_3d( x(m:m+1), y(n:n+1), z(l:l+1), w(m:m+1,n:n+1,l:l+1), inter_p, inter_v(3) )  

        k2=dt*inter_v(1)
        l2=dt*inter_v(2)
        n2=dt*inter_v(3)
        x1=trajx(i)+0.5*(k1+k2)
        y1=trajy(i)+0.5*(l1+l2)
        z1=trajz(i)+0.5*(n1+n2)

!-- 計算した traj が領域内に存在しているか確認
        if(x1<x(1).or.x1>x(nx).or.y1<y(1).or.y1>y(ny).or.z1<z(1).or.z1>z(nz))then
           write(*,*) "****** ERROR ******"
           write(*,*) "This point does not exist in the region."
           write(*,*) "Stop.!!"
!-- 計算領域外に出た場合は, 以降のステップにすべて defun を代入
           do id=i,step-1
              trajx(id+1)=defun
              trajy(id+1)=defun
              trajz(id+1)=defun
           end do
           exit
        end if
        trajx(i+1)=x1
        trajy(i+1)=y1
        trajz(i+1)=z1
     end do


  case ('ME1')
     do i=1,step-1
        inter_p=(/trajx(i), trajy(i), trajz(i)/)
        call interpo_search_3d( x, y, z, trajx(i), trajy(i), trajz(i), m, n, l )  
        call interpolation_3d( x(m:m+1), y(n:n+1), z(l:l+1), u(m:m+1,n:n+1,l:l+1), inter_p, inter_v(1) )  
        call interpolation_3d( x(m:m+1), y(n:n+1), z(l:l+1), v(m:m+1,n:n+1,l:l+1), inter_p, inter_v(2) )  
        call interpolation_3d( x(m:m+1), y(n:n+1), z(l:l+1), w(m:m+1,n:n+1,l:l+1), inter_p, inter_v(3) )  

        k1=dt*inter_v(1)
        l1=dt*inter_v(2)
        n1=dt*inter_v(3)
!-- 一時的な流跡線の位置を計算
        x1=trajx(i)+0.5*k1
        y1=trajy(i)+0.5*l1
        z1=trajz(i)+0.5*n1
!-- この点での速度場を内挿
        inter_p=(/x1, y1, z1/)
        call interpo_search_3d( x, y, z, x1, y1, z1, m, n, l )  
        call interpolation_3d( x(m:m+1), y(n:n+1), z(l:l+1), u(m:m+1,n:n+1,l:l+1), inter_p, inter_v(1) )  
        call interpolation_3d( x(m:m+1), y(n:n+1), z(l:l+1), v(m:m+1,n:n+1,l:l+1), inter_p, inter_v(2) )  
        call interpolation_3d( x(m:m+1), y(n:n+1), z(l:l+1), w(m:m+1,n:n+1,l:l+1), inter_p, inter_v(3) )  

        k2=dt*inter_v(1)
        l2=dt*inter_v(2)
        n2=dt*inter_v(3)
        x1=trajx(i)+k2
        y1=trajy(i)+l2
        z1=trajz(i)+n2

!-- 計算した traj が領域内に存在しているか確認
        if(x1<x(1).or.x1>x(nx).or.y1<y(1).or.y1>y(ny).or.z1<z(1).or.z1>z(nz))then
           write(*,*) "****** ERROR ******"
           write(*,*) "This point does not exist in the region."
           write(*,*) "Stop.!!"
!-- 計算領域外に出た場合は, 以降のステップにすべて defun を代入
           do id=i,step-1
              trajx(id+1)=defun
              trajy(id+1)=defun
              trajz(id+1)=defun
           end do
           exit
        end if
        trajx(i+1)=x1
        trajy(i+1)=y1
        trajz(i+1)=z1
     end do



  case ('RK3')
     do i=1,step-1
        inter_p=(/trajx(i), trajy(i), trajz(i)/)
        call interpo_search_3d( x, y, z, trajx(i), trajy(i), trajz(i), m, n, l )  
        call interpolation_3d( x(m:m+1), y(n:n+1), z(l:l+1), u(m:m+1,n:n+1,l:l+1), inter_p, inter_v(1) )  
        call interpolation_3d( x(m:m+1), y(n:n+1), z(l:l+1), v(m:m+1,n:n+1,l:l+1), inter_p, inter_v(2) )  
        call interpolation_3d( x(m:m+1), y(n:n+1), z(l:l+1), w(m:m+1,n:n+1,l:l+1), inter_p, inter_v(3) )  

        k1=dt*inter_v(1)
        l1=dt*inter_v(2)
        n1=dt*inter_v(3)
!-- 一時的な流跡線の位置を計算
        x1=trajx(i)+0.5*k1
        y1=trajy(i)+0.5*l1
        z1=trajz(i)+0.5*n1
!-- この点での速度場を内挿
        inter_p=(/x1, y1, z1/)
        call interpo_search_3d( x, y, z, x1, y1, z1, m, n, l )  
        call interpolation_3d( x(m:m+1), y(n:n+1), z(l:l+1), u(m:m+1,n:n+1,l:l+1), inter_p, inter_v(1) )  
        call interpolation_3d( x(m:m+1), y(n:n+1), z(l:l+1), v(m:m+1,n:n+1,l:l+1), inter_p, inter_v(2) )  
        call interpolation_3d( x(m:m+1), y(n:n+1), z(l:l+1), w(m:m+1,n:n+1,l:l+1), inter_p, inter_v(3) )  

        k2=dt*inter_v(1)
        l2=dt*inter_v(2)
        n2=dt*inter_v(3)
        x1=trajx(i)-k1+2.0*k2
        y1=trajy(i)-l1+2.0*l2
        z1=trajz(i)-n1+2.0*n2

!-- この点での速度場を内挿
        inter_p=(/x1, y1, z1/)
        call interpo_search_3d( x, y, z, x1, y1, z1, m, n, l )  
        call interpolation_3d( x(m:m+1), y(n:n+1), z(l:l+1), u(m:m+1,n:n+1,l:l+1), inter_p, inter_v(1) )  
        call interpolation_3d( x(m:m+1), y(n:n+1), z(l:l+1), v(m:m+1,n:n+1,l:l+1), inter_p, inter_v(2) )  
        call interpolation_3d( x(m:m+1), y(n:n+1), z(l:l+1), w(m:m+1,n:n+1,l:l+1), inter_p, inter_v(3) )  

        k3=dt*inter_v(1)
        l3=dt*inter_v(2)
        n3=dt*inter_v(3)
        x1=trajx(i)+(1.0/6.0)*(k1+4.0*k2+k3)
        y1=trajy(i)+(1.0/6.0)*(l1+4.0*l2+l3)
        z1=trajz(i)+(1.0/6.0)*(n1+4.0*n2+n3)

!-- 計算した traj が領域内に存在しているか確認
        if(x1<x(1).or.x1>x(nx).or.y1<y(1).or.y1>y(ny).or.z1<z(1).or.z1>z(nz))then
           write(*,*) "****** ERROR ******"
           write(*,*) "This point does not exist in the region."
           write(*,*) "Stop.!!"
!-- 計算領域外に出た場合は, 以降のステップにすべて defun を代入
           do id=i,step-1
              trajx(id+1)=defun
              trajy(id+1)=defun
              trajz(id+1)=defun
           end do
           exit
        end if
        trajx(i+1)=x1
        trajy(i+1)=y1
        trajz(i+1)=z1
     end do


  case ('RK4')
     do i=1,step-1
        inter_p=(/trajx(i), trajy(i), trajz(i)/)
        call interpo_search_3d( x, y, z, trajx(i), trajy(i), trajz(i), m, n, l )  
        call interpolation_3d( x(m:m+1), y(n:n+1), z(l:l+1), u(m:m+1,n:n+1,l:l+1), inter_p, inter_v(1) )  
        call interpolation_3d( x(m:m+1), y(n:n+1), z(l:l+1), v(m:m+1,n:n+1,l:l+1), inter_p, inter_v(2) )  
        call interpolation_3d( x(m:m+1), y(n:n+1), z(l:l+1), w(m:m+1,n:n+1,l:l+1), inter_p, inter_v(3) )  

        k1=dt*inter_v(1)
        l1=dt*inter_v(2)
        n1=dt*inter_v(3)
!-- 一時的な流跡線の位置を計算
        x1=trajx(i)+0.5*k1
        y1=trajy(i)+0.5*l1
        z1=trajz(i)+0.5*n1
!-- この点での速度場を内挿
        inter_p=(/x1, y1, z1/)
        call interpo_search_3d( x, y, z, x1, y1, z1, m, n, l )  
        call interpolation_3d( x(m:m+1), y(n:n+1), z(l:l+1), u(m:m+1,n:n+1,l:l+1), inter_p, inter_v(1) )  
        call interpolation_3d( x(m:m+1), y(n:n+1), z(l:l+1), v(m:m+1,n:n+1,l:l+1), inter_p, inter_v(2) )  
        call interpolation_3d( x(m:m+1), y(n:n+1), z(l:l+1), w(m:m+1,n:n+1,l:l+1), inter_p, inter_v(3) )  

        k2=dt*inter_v(1)
        l2=dt*inter_v(2)
        n2=dt*inter_v(3)
        x1=trajx(i)+0.5*k2
        y1=trajy(i)+0.5*l2
        z1=trajz(i)+0.5*n2

!-- この点での速度場を内挿
        inter_p=(/x1, y1, z1/)
        call interpo_search_3d( x, y, z, x1, y1, z1, m, n, l )  
        call interpolation_3d( x(m:m+1), y(n:n+1), z(l:l+1), u(m:m+1,n:n+1,l:l+1), inter_p, inter_v(1) )  
        call interpolation_3d( x(m:m+1), y(n:n+1), z(l:l+1), v(m:m+1,n:n+1,l:l+1), inter_p, inter_v(2) )  
        call interpolation_3d( x(m:m+1), y(n:n+1), z(l:l+1), w(m:m+1,n:n+1,l:l+1), inter_p, inter_v(3) )  

        k3=dt*inter_v(1)
        l3=dt*inter_v(2)
        n3=dt*inter_v(3)
        x1=trajx(i)+k3
        y1=trajy(i)+l3
        z1=trajz(i)+n3

!-- この点での速度場を内挿
        inter_p=(/x1, y1, z1/)
        call interpo_search_3d( x, y, z, x1, y1, z1, m, n, l )  
        call interpolation_3d( x(m:m+1), y(n:n+1), z(l:l+1), u(m:m+1,n:n+1,l:l+1), inter_p, inter_v(1) )  
        call interpolation_3d( x(m:m+1), y(n:n+1), z(l:l+1), v(m:m+1,n:n+1,l:l+1), inter_p, inter_v(2) )  
        call interpolation_3d( x(m:m+1), y(n:n+1), z(l:l+1), w(m:m+1,n:n+1,l:l+1), inter_p, inter_v(3) )  

        k4=dt*inter_v(1)
        l4=dt*inter_v(2)
        n4=dt*inter_v(3)
        x1=trajx(i)+(1.0/6.0)*(k1+2.0*k2+2.0*k3+k4)
        y1=trajy(i)+(1.0/6.0)*(l1+2.0*l2+2.0*l3+l4)
        z1=trajz(i)+(1.0/6.0)*(n1+2.0*n2+2.0*n3+n4)

!-- 計算した traj が領域内に存在しているか確認
        if(x1<x(1).or.x1>x(nx).or.y1<y(1).or.y1>y(ny).or.z1<z(1).or.z1>z(nz))then
           write(*,*) "****** ERROR ******"
           write(*,*) "This point does not exist in the region."
           write(*,*) "Stop.!!"
!-- 計算領域外に出た場合は, 以降のステップにすべて defun を代入
           do id=i,step-1
              trajx(id+1)=defun
              trajy(id+1)=defun
              trajz(id+1)=defun
           end do
           exit
        end if
        trajx(i+1)=x1
        trajy(i+1)=y1
        trajz(i+1)=z1
     end do



  end select

end subroutine Forward_Traject_3d