| Class | Trajectory |
| In: |
trajectory.f90
|
トラジェクトリ解析を行うルーチン
| Subroutine : | |||
| dt : | real, intent(in)
| ||
| stime : | real, intent(in)
| ||
| step : | integer, intent(in)
| ||
| ini_x(:) : | real, intent(in)
| ||
| ini_y(size(ini_x)) : | real, intent(in)
| ||
| t(:) : | real, intent(in)
| ||
| x(:) : | real, intent(in)
| ||
| y(:) : | real, intent(in)
| ||
| u(size(x),size(y),size(t)) : | real, intent(in)
| ||
| v(size(x),size(y),size(t)) : | real, intent(in)
| ||
| trajx(step,size(ini_x)) : | real, intent(inout)
| ||
| trajy(step,size(ini_x)) : | real, intent(inout)
| ||
| FTF(size(ini_x)) : | logical, intent(inout)
| ||
| opt : | character(*), intent(in), optional
| ||
| undef : | real, intent(in), optional
| ||
| stdopt : | logical, intent(in), optional
|
オフラインでの後方流跡線解析ルーチン. Stream_Line_2d を毎ステップ呼び出すが, その際に代入される流速場は, 引数として与えられる風速場を時間方向に線形内挿したもので毎時 u,v が更新される.
subroutine Backward_Traject_2d( dt, stime, step, ini_x, ini_y, t, x, y, u, v, trajx, trajy, FTF, opt, undef, stdopt )
! オフラインでの後方流跡線解析ルーチン.
! Stream_Line_2d を毎ステップ呼び出すが, その際に代入される流速場は,
! 引数として与えられる風速場を時間方向に線形内挿したもので毎時 u,v が更新される.
implicit none
real, intent(in) :: dt ! 計算する時間間隔 (正値) [s]
real, intent(in) :: stime ! 初期位置での時刻 [s]
integer, intent(in) :: step ! 計算を行うステップ数
real, intent(in) :: ini_x(:) ! 初期位置 x [m] 複数個設定可
real, intent(in) :: ini_y(size(ini_x)) ! 初期位置 y [m] 複数個設定可
real, intent(in) :: t(:) ! 速度場データのある時刻 [s]
! 時刻は t(1) > stime, t(size(t)) < stime - dt*step
real, intent(in) :: x(:) ! x 方向の座標値 [m]
real, intent(in) :: y(:) ! y 方向の座標値 [m]
real, intent(in) :: u(size(x),size(y),size(t)) ! x 方向の風速成分 [m/s]
real, intent(in) :: v(size(x),size(y),size(t)) ! y 方向の風速成分 [m/s]
real, intent(inout) :: trajx(step,size(ini_x)) ! トラジェクトリの位置座標 x [m]
real, intent(inout) :: trajy(step,size(ini_x)) ! トラジェクトリの位置座標 y [m]
logical, intent(inout) :: FTF(size(ini_x)) ! 各パーセルが領域外に出たか.
! 出て入いれば true,
! 出ていなければ false.
character(*), intent(in), optional :: opt ! 時間積分スキーム
real, intent(in), optional :: undef ! 未定義
logical, intent(in), optional :: stdopt ! 計算領域外警告除外フラグ
! デフォルトは false (表示する)
character(3) :: option
integer :: nx, ny, nt, nc, i, j, k, i_trajt, ri_trajt
integer :: defun_i
real :: inter_u(size(x),size(y)), inter_v(size(x),size(y))
real :: trajt(step), r_t(size(t))
real :: defun
logical :: stderr
nx=size(x)
ny=size(y)
nt=size(t)
nc=size(ini_x)
! FTF=.false.
trajt=(/((stime-dt*(i-1)),i=1,step)/) ! 後方なので -dt で計算.
do i=1,nt
r_t(i)=t(nt-i+1)
end do
do i=1,nc
trajx(1,i)=ini_x(i)
trajy(1,i)=ini_y(i)
end do
!-- 速度場データが stime 以前にない場合, 内挿できないので, エラー
if(trajt(1)>t(1))then
write(*,*) "#### ERROR ####"
write(*,*) "stime must be less than t(1)."
write(*,*) "STOP!!"
stop
end if
if(trajt(step)<t(nt))then
write(*,*) "#### ERROR ####"
write(*,*) "trajt(step) must be more than t(nt)."
write(*,*) "STOP!!"
stop
end if
if(present(opt))then
option=opt
else
option='EU1'
end if
if(present(undef))then
defun=undef
defun_i=int(undef)
else
defun_i=-2147483648
defun=real(defun_i)
end if
do j=1,nc
do i=2,step
trajx(i,j)=defun
trajy(i,j)=defun
end do
end do
if(present(stdopt))then
stderr=stdopt
else
stderr=.false.
end if
do i=1,step-1
call interpo_search_1d( r_t, trajt(i), ri_trajt )
! 流跡線時間は trajt 基準でスタート.
! --->--->--->--->---> 時間の流れ --->--->--->--->--->--->
! ------ r_t(ri_trajt) ----------- r_t(ri_trajt+1) -------
! | |
! ---- t(nt-ri_trajt+1) ----------- t(nt-ri_trajt) -------
! ---------------------- trajt(i) ------------------------
! 速度変数は全て逆向き時間で入っている.
! forward で内挿する時間領域が i_trajt, i_trajt+1 の場合,
! backward では i_trajt, i_trajt-1 となる.
! i_trajt = nt - ri_trajt + 1
i_trajt=nt-ri_trajt+1
!-- 速度場の内挿
if(t(i_trajt)==trajt(i))then
if(present(undef))then
do k=1,ny
do j=1,nx
if(u(j,k,i_trajt)/=undef.and.v(j,k,i_trajt)/=undef)then
inter_u(j,k)=u(j,k,i_trajt)
inter_v(j,k)=v(j,k,i_trajt)
else
inter_u(j,k)=undef
inter_v(j,k)=undef
end if
end do
end do
else
do k=1,ny
do j=1,nx
inter_u(j,k)=u(j,k,i_trajt)
inter_v(j,k)=v(j,k,i_trajt)
end do
end do
end if
else
if(present(undef))then
!$omp parallel default(shared)
!$omp do schedule(dynamic) private(j,k)
do k=1,ny
do j=1,nx
if(u(j,k,i_trajt)/=undef.and.u(j,k,i_trajt-1)/=undef.and. v(j,k,i_trajt)/=undef.and.v(j,k,i_trajt-1)/=undef)then
call interpolation_1d( (/t(i_trajt), t(i_trajt-1)/), (/u(j,k,i_trajt), u(j,k,i_trajt-1)/), trajt(i), inter_u(j,k) )
call interpolation_1d( (/t(i_trajt), t(i_trajt-1)/), (/v(j,k,i_trajt), v(j,k,i_trajt-1)/), trajt(i), inter_v(j,k) )
else
inter_u(j,k)=undef
inter_v(j,k)=undef
end if
end do
end do
!$omp end do
!$omp end parallel
else
!$omp parallel default(shared)
!$omp do schedule(runtime) private(j,k)
do k=1,ny
do j=1,nx
call interpolation_1d( (/t(i_trajt), t(i_trajt-1)/), (/u(j,k,i_trajt), u(j,k,i_trajt-1)/), trajt(i), inter_u(j,k) )
call interpolation_1d( (/t(i_trajt), t(i_trajt-1)/), (/v(j,k,i_trajt), v(j,k,i_trajt-1)/), trajt(i), inter_v(j,k) )
end do
end do
!$omp end do
!$omp end parallel
end if
end if
if(present(undef))then
!$omp parallel default(shared)
!$omp do schedule(dynamic) private(k)
do k=1,nc
if(FTF(k).eqv..false.)then
call Stream_Line_2d( -dt, 2, trajx(i,k), trajy(i,k), x, y, inter_u, inter_v, trajx(i:i+1,k), trajy(i:i+1,k), FTF(k), opt=option, undef=undef, stdopt=stderr )
end if
end do
!$omp end do
!$omp end parallel
else
!$omp parallel default(shared)
!$omp do schedule(dynamic) private(k)
do k=1,nc
if(FTF(k).eqv..false.)then
call Stream_Line_2d( -dt, 2, trajx(i,k), trajy(i,k), x, y, inter_u, inter_v, trajx(i:i+1,k), trajy(i:i+1,k), FTF(k), opt=option, stdopt=stderr )
end if
end do
!$omp end do
!$omp end parallel
end if
end do
end subroutine Backward_Traject_2d
| Subroutine : | |||
| dt : | real, intent(in)
| ||
| stime : | real, intent(in)
| ||
| step : | integer, intent(in)
| ||
| ini_x(:) : | real, intent(in)
| ||
| ini_y(size(ini_x)) : | real, intent(in)
| ||
| ini_z(size(ini_x)) : | real, intent(in)
| ||
| t(:) : | real, intent(in)
| ||
| x(:) : | real, intent(in)
| ||
| y(:) : | real, intent(in)
| ||
| z(:) : | real, intent(in)
| ||
| u(size(x),size(y),size(z),size(t)) : | real, intent(in)
| ||
| v(size(x),size(y),size(z),size(t)) : | real, intent(in)
| ||
| w(size(x),size(y),size(z),size(t)) : | real, intent(in)
| ||
| trajx(step,(size(ini_x))) : | real, intent(inout)
| ||
| trajy(step,(size(ini_x))) : | real, intent(inout)
| ||
| trajz(step,(size(ini_x))) : | real, intent(inout)
| ||
| FTF(size(ini_x)) : | logical, intent(inout)
| ||
| opt : | character(*), intent(in), optional
| ||
| undef : | real, intent(in), optional
| ||
| stdopt : | logical, intent(in), optional
|
オフラインでの後方流跡線解析ルーチン. Stream_Line_3d を毎ステップ呼び出すが, その際に代入される流速場は, 引数として与えられる風速場を時間方向に線形内挿したもので毎時 u,v が更新される.
subroutine Backward_Traject_3d( dt, stime, step, ini_x, ini_y, ini_z, t, x, y, z, u, v, w, trajx, trajy, trajz, FTF, opt, undef, stdopt )
! オフラインでの後方流跡線解析ルーチン.
! Stream_Line_3d を毎ステップ呼び出すが, その際に代入される流速場は,
! 引数として与えられる風速場を時間方向に線形内挿したもので毎時 u,v が更新される.
implicit none
real, intent(in) :: dt ! 計算する時間間隔 [s]
real, intent(in) :: stime ! 初期位置での時刻 [s]
integer, intent(in) :: step ! 計算を行うステップ数
real, intent(in) :: ini_x(:) ! 初期位置 x [m] 複数個設定可
real, intent(in) :: ini_y(size(ini_x)) ! 初期位置 y [m] 複数個設定可
real, intent(in) :: ini_z(size(ini_x)) ! 初期位置 z [m] 複数個設定可
real, intent(in) :: t(:) ! 速度場データのある時刻 [s]
! 時刻は t(size(t)) > stime, t(1) < stime - dt*step
real, intent(in) :: x(:) ! x 方向の座標値 [m]
real, intent(in) :: y(:) ! y 方向の座標値 [m]
real, intent(in) :: z(:) ! z 方向の座標値 [m]
real, intent(in) :: u(size(x),size(y),size(z),size(t)) ! x 方向の風速成分 [m/s]
real, intent(in) :: v(size(x),size(y),size(z),size(t)) ! y 方向の風速成分 [m/s]
real, intent(in) :: w(size(x),size(y),size(z),size(t)) ! z 方向の風速成分 [m/s]
real, intent(inout) :: trajx(step,(size(ini_x))) ! トラジェクトリの位置座標 x [m]
real, intent(inout) :: trajy(step,(size(ini_x))) ! トラジェクトリの位置座標 y [m]
real, intent(inout) :: trajz(step,(size(ini_x))) ! トラジェクトリの位置座標 z [m]
logical, intent(inout) :: FTF(size(ini_x)) ! パーセルが領域外に出たか.
! 出たら true,
! 出なければ false.
character(*), intent(in), optional :: opt ! 時間積分スキーム
real, intent(in), optional :: undef ! 未定義
logical, intent(in), optional :: stdopt ! 計算領域外警告除外フラグ
! デフォルトは false (表示する)
character(3) :: option
integer :: nx, ny, nz, nt, nc, i, j, k, l, i_trajt, ri_trajt
integer :: defun_i
real :: inter_u(size(x),size(y),size(z))
real :: inter_v(size(x),size(y),size(z))
real :: inter_w(size(x),size(y),size(z))
real :: trajt(step), r_t(size(t))
real :: defun
logical :: stderr
nx=size(x)
ny=size(y)
nz=size(z)
nt=size(t)
nc=size(ini_x)
! FTF=.false.
trajt=(/((stime-dt*(i-1)),i=1,step)/) ! 後方なので -dt で計算.
do i=1,nt
r_t(i)=t(nt-i+1)
end do
do i=1,nc
trajx(1,i)=ini_x(i)
trajy(1,i)=ini_y(i)
trajz(1,i)=ini_z(i)
end do
!-- 速度場データが stime 以前にない場合, 内挿できないので, エラー
if(trajt(1)>t(1))then
write(*,*) "#### ERROR ####"
write(*,*) "stime must be less than t(1)."
write(*,*) "STOP!!"
stop
end if
if(trajt(step)<t(nt))then
write(*,*) "#### ERROR ####"
write(*,*) "trajt(step) must be more than t(nt)."
write(*,*) "STOP!!"
stop
end if
if(present(opt))then
option=opt
else
option='EU1'
end if
if(present(undef))then
defun=undef
defun_i=int(undef)
else
defun_i=-2147483648
defun=real(defun_i)
end if
do j=1,nc
do i=2,step
trajx(i,j)=defun
trajy(i,j)=defun
trajz(i,j)=defun
end do
end do
if(present(stdopt))then
stderr=stdopt
else
stderr=.false.
end if
!write(*,*) "step check", step
do i=1,step-1
call interpo_search_1d( r_t, trajt(i), ri_trajt )
! --->--->--->--->---> 時間の流れ --->--->--->--->--->--->
! ------ r_t(ri_trajt) ----------- r_t(ri_trajt+1) -------
! | |
! ---- t(nt-ri_trajt+1) ----------- t(nt-ri_trajt) -------
! ---------------------- trajt(i) ------------------------
! 速度変数は全て逆向き時間で入っている.
! forward で内挿する時間領域が i_trajt, i_trajt+1 の場合,
! backward では i_trajt, i_trajt-1 となる.
! i_trajt = nt - ri_trajt + 1
i_trajt=nt-ri_trajt+1
!write(*,*) "step_traj", trajt(i)
!-- 速度場の内挿
if(t(i_trajt)==trajt(i))then
if(present(undef))then
do l=1,nz
do k=1,ny
do j=1,nx
if(u(j,k,l,i_trajt)/=undef.and.v(j,k,l,i_trajt)/=undef .and.w(j,k,l,i_trajt)/=undef)then
inter_u(j,k,l)=u(j,k,l,i_trajt)
inter_v(j,k,l)=v(j,k,l,i_trajt)
inter_w(j,k,l)=w(j,k,l,i_trajt)
else
inter_u(j,k,l)=undef
inter_v(j,k,l)=undef
inter_w(j,k,l)=undef
end if
end do
end do
end do
else
do l=1,nz
do k=1,ny
do j=1,nx
inter_u(j,k,l)=u(j,k,l,i_trajt)
inter_v(j,k,l)=v(j,k,l,i_trajt)
inter_w(j,k,l)=w(j,k,l,i_trajt)
end do
end do
end do
end if
else
!write(*,*) "interpo_check, trajt", t, trajt(i), i_trajt
if(present(undef))then
!$omp parallel default(shared)
!$omp do schedule(dynamic) private(j,k,l)
do l=1,nz
do k=1,ny
do j=1,nx
if(u(j,k,l,i_trajt)/=undef.and.u(j,k,l,i_trajt-1)/=undef.and. v(j,k,l,i_trajt)/=undef.and.v(j,k,l,i_trajt-1)/=undef.and. w(j,k,l,i_trajt)/=undef.and.w(j,k,l,i_trajt-1)/=undef)then
call interpolation_1d( (/t(i_trajt), t(i_trajt-1)/), (/u(j,k,l,i_trajt), u(j,k,l,i_trajt-1)/), trajt(i), inter_u(j,k,l) )
call interpolation_1d( (/t(i_trajt), t(i_trajt-1)/), (/v(j,k,l,i_trajt), v(j,k,l,i_trajt-1)/), trajt(i), inter_v(j,k,l) )
call interpolation_1d( (/t(i_trajt), t(i_trajt-1)/), (/w(j,k,l,i_trajt), w(j,k,l,i_trajt-1)/), trajt(i), inter_w(j,k,l) )
else
inter_u(j,k,l)=undef
inter_v(j,k,l)=undef
inter_w(j,k,l)=undef
end if
end do
end do
end do
!$omp end do
!$omp end parallel
else
!$omp parallel default(shared)
!$omp do schedule(runtime) private(j,k,l)
do l=1,nz
do k=1,ny
do j=1,nx
call interpolation_1d( (/t(i_trajt), t(i_trajt-1)/), (/u(j,k,l,i_trajt), u(j,k,l,i_trajt-1)/), trajt(i), inter_u(j,k,l) )
call interpolation_1d( (/t(i_trajt), t(i_trajt-1)/), (/v(j,k,l,i_trajt), v(j,k,l,i_trajt-1)/), trajt(i), inter_v(j,k,l) )
call interpolation_1d( (/t(i_trajt), t(i_trajt-1)/), (/w(j,k,l,i_trajt), w(j,k,l,i_trajt-1)/), trajt(i), inter_w(j,k,l) )
end do
end do
end do
!$omp end do
!$omp end parallel
end if
end if
if(present(undef))then
!$omp parallel default(shared)
!$omp do schedule(dynamic) private(k)
do k=1,nc
if(FTF(k).eqv..false.)then
call Stream_Line_3d( -dt, 2, trajx(i,k), trajy(i,k), trajz(i,k), x, y, z, inter_u, inter_v, inter_w, trajx(i:i+1,k), trajy(i:i+1,k), trajz(i:i+1,k), FTF(k), opt=option, undef=undef, stdopt=stderr )
end if
end do
!$omp end do
!$omp end parallel
else
!$omp parallel default(shared)
!$omp do schedule(dynamic) private(k)
do k=1,nc
if(FTF(k).eqv..false.)then
call Stream_Line_3d( -dt, 2, trajx(i,k), trajy(i,k), trajz(i,k), x, y, z, inter_u, inter_v, inter_w, trajx(i:i+1,k), trajy(i:i+1,k), trajz(i:i+1,k), FTF(k), opt=option, stdopt=stderr )
end if
end do
!$omp end do
!$omp end parallel
end if
end do
end subroutine Backward_Traject_3d
| Subroutine : | |||
| dt : | real, intent(in)
| ||
| stime : | real, intent(in)
| ||
| step : | integer, intent(in)
| ||
| ini_x(:) : | real, intent(in)
| ||
| ini_y(size(ini_x)) : | real, intent(in)
| ||
| t(:) : | real, intent(in)
| ||
| x(:) : | real, intent(in)
| ||
| y(:) : | real, intent(in)
| ||
| u(size(x),size(y),size(t)) : | real, intent(in)
| ||
| v(size(x),size(y),size(t)) : | real, intent(in)
| ||
| trajx(step,size(ini_x)) : | real, intent(inout)
| ||
| trajy(step,size(ini_x)) : | real, intent(inout)
| ||
| FTF(size(ini_x)) : | logical, intent(inout)
| ||
| opt : | character(*), intent(in), optional
| ||
| undef : | real, intent(in), optional
| ||
| stdopt : | logical, intent(in), optional
|
オフラインでの流跡線解析ルーチン. Stream_Line_2d を毎ステップ呼び出すが, その際に代入される流速場は, 引数として与えられる風速場を時間方向に線形内挿したもので毎時 u,v が更新される.
subroutine Forward_Traject_2d( dt, stime, step, ini_x, ini_y, t, x, y, u, v, trajx, trajy, FTF, opt, undef, stdopt )
! オフラインでの流跡線解析ルーチン.
! Stream_Line_2d を毎ステップ呼び出すが, その際に代入される流速場は,
! 引数として与えられる風速場を時間方向に線形内挿したもので毎時 u,v が更新される.
implicit none
real, intent(in) :: dt ! 計算する時間間隔 [s]
real, intent(in) :: stime ! 初期位置での時刻 [s]
integer, intent(in) :: step ! 計算を行うステップ数
real, intent(in) :: ini_x(:) ! 初期位置 x [m] 複数個設定可
real, intent(in) :: ini_y(size(ini_x)) ! 初期位置 y [m] 複数個設定可
real, intent(in) :: t(:) ! 速度場データのある時刻 [s]
real, intent(in) :: x(:) ! x 方向の座標値 [m]
real, intent(in) :: y(:) ! y 方向の座標値 [m]
real, intent(in) :: u(size(x),size(y),size(t)) ! x 方向の風速成分 [m/s]
real, intent(in) :: v(size(x),size(y),size(t)) ! y 方向の風速成分 [m/s]
real, intent(inout) :: trajx(step,size(ini_x)) ! トラジェクトリの位置座標 x [m]
real, intent(inout) :: trajy(step,size(ini_x)) ! トラジェクトリの位置座標 y [m]
logical, intent(inout) :: FTF(size(ini_x)) ! 各パーセルが領域外に出たか.
! 出て入いれば true,
! 出ていなければ false
character(*), intent(in), optional :: opt ! 時間積分スキーム
real, intent(in), optional :: undef ! 未定義
logical, intent(in), optional :: stdopt ! 計算領域外警告除外フラグ
! デフォルトは false (表示する)
character(3) :: option
integer :: nx, ny, nt, nc, i, j, k, i_trajt
integer :: defun_i
real :: inter_u(size(x),size(y)), inter_v(size(x),size(y))
real :: trajt(step)
real :: defun
logical :: stderr
nx=size(x)
ny=size(y)
nt=size(t)
nc=size(ini_x)
! FTF=.false.
trajt=(/((stime+dt*(i-1)),i=1,step)/)
do i=1,nc
trajx(1,i)=ini_x(i)
trajy(1,i)=ini_y(i)
end do
!-- 速度場データが stime 以前にない場合, 内挿できないので, エラー
if(trajt(1)<t(1))then
write(*,*) "#### ERROR ####"
write(*,*) "stime must be more than t(1)."
write(*,*) "STOP!!"
stop
end if
if(trajt(step)>t(nt))then
write(*,*) "#### ERROR ####"
write(*,*) "trajt(step) must be less than t(nt)."
write(*,*) "STOP!!"
stop
end if
if(present(opt))then
option=opt
else
option='EU1'
end if
if(present(undef))then
defun=undef
defun_i=int(undef)
else
defun_i=-2147483648
defun=real(defun_i)
end if
do j=1,nc
do i=2,step
trajx(i,j)=defun
trajy(i,j)=defun
end do
end do
if(present(stdopt))then
stderr=stdopt
else
stderr=.false.
end if
do i=1,step-1
call interpo_search_1d( t, trajt(i), i_trajt )
!-- 速度場の内挿
if(t(i_trajt)==trajt(i))then
if(present(undef))then
do k=1,ny
do j=1,nx
if(u(j,k,i_trajt)/=undef.and.v(j,k,i_trajt)/=undef)then
inter_u(j,k)=u(j,k,i_trajt)
inter_v(j,k)=v(j,k,i_trajt)
else
inter_u(j,k)=undef
inter_v(j,k)=undef
end if
end do
end do
else
do k=1,ny
do j=1,nx
inter_u(j,k)=u(j,k,i_trajt)
inter_v(j,k)=v(j,k,i_trajt)
end do
end do
end if
else
if(present(undef))then
!$omp parallel default(shared)
!$omp do schedule(dynamic) private(j,k)
do k=1,ny
do j=1,nx
if(u(j,k,i_trajt)/=undef.and.u(j,k,i_trajt+1)/=undef.and. v(j,k,i_trajt)/=undef.and.v(j,k,i_trajt+1)/=undef)then
call interpolation_1d( (/t(i_trajt), t(i_trajt+1)/), (/u(j,k,i_trajt), u(j,k,i_trajt+1)/), trajt(i), inter_u(j,k) )
call interpolation_1d( (/t(i_trajt), t(i_trajt+1)/), (/v(j,k,i_trajt), v(j,k,i_trajt+1)/), trajt(i), inter_v(j,k) )
else
inter_u(j,k)=undef
inter_v(j,k)=undef
end if
end do
end do
!$omp end do
!$omp end parallel
else
!$omp parallel default(shared)
!$omp do schedule(runtime) private(j,k)
do k=1,ny
do j=1,nx
call interpolation_1d( (/t(i_trajt), t(i_trajt+1)/), (/u(j,k,i_trajt), u(j,k,i_trajt+1)/), trajt(i), inter_u(j,k) )
call interpolation_1d( (/t(i_trajt), t(i_trajt+1)/), (/v(j,k,i_trajt), v(j,k,i_trajt+1)/), trajt(i), inter_v(j,k) )
end do
end do
!$omp end do
!$omp end parallel
end if
end if
if(present(undef))then
!$omp parallel default(shared)
!$omp do schedule(dynamic) private(k)
do k=1,nc
if(FTF(k).eqv..false.)then
call Stream_Line_2d( dt, 2, trajx(i,k), trajy(i,k), x, y, inter_u, inter_v, trajx(i:i+1,k), trajy(i:i+1,k), FTF(k), opt=option, undef=undef, stdopt=stderr )
end if
end do
!$omp end do
!$omp end parallel
else
!$omp parallel default(shared)
!$omp do schedule(dynamic) private(k)
do k=1,nc
if(FTF(k).eqv..false.)then
call Stream_Line_2d( dt, 2, trajx(i,k), trajy(i,k), x, y, inter_u, inter_v, trajx(i:i+1,k), trajy(i:i+1,k), FTF(k), opt=option, stdopt=stderr )
end if
end do
!$omp end do
!$omp end parallel
end if
end do
end subroutine Forward_Traject_2d
| Subroutine : | |||
| dt : | real, intent(in)
| ||
| stime : | real, intent(in)
| ||
| step : | integer, intent(in)
| ||
| ini_x(:) : | real, intent(in)
| ||
| ini_y(size(ini_x)) : | real, intent(in)
| ||
| ini_z(size(ini_x)) : | real, intent(in)
| ||
| t(:) : | real, intent(in)
| ||
| x(:) : | real, intent(in)
| ||
| y(:) : | real, intent(in)
| ||
| z(:) : | real, intent(in)
| ||
| u(size(x),size(y),size(z),size(t)) : | real, intent(in)
| ||
| v(size(x),size(y),size(z),size(t)) : | real, intent(in)
| ||
| w(size(x),size(y),size(z),size(t)) : | real, intent(in)
| ||
| trajx(step,(size(ini_x))) : | real, intent(inout)
| ||
| trajy(step,(size(ini_x))) : | real, intent(inout)
| ||
| trajz(step,(size(ini_x))) : | real, intent(inout)
| ||
| FTF(size(ini_x)) : | logical, intent(inout)
| ||
| opt : | character(*), intent(in), optional
| ||
| undef : | real, intent(in), optional
| ||
| stdopt : | logical, intent(in), optional
|
オフラインでの前方流跡線解析ルーチン. Stream_Line_3d を毎ステップ呼び出すが, その際に代入される流速場は, 引数として与えられる風速場を時間方向に線形内挿したもので毎時 u,v が更新される.
subroutine Forward_Traject_3d( dt, stime, step, ini_x, ini_y, ini_z, t, x, y, z, u, v, w, trajx, trajy, trajz, FTF, opt, undef, stdopt )
! オフラインでの前方流跡線解析ルーチン.
! Stream_Line_3d を毎ステップ呼び出すが, その際に代入される流速場は,
! 引数として与えられる風速場を時間方向に線形内挿したもので毎時 u,v が更新される.
implicit none
real, intent(in) :: dt ! 計算する時間間隔 [s]
real, intent(in) :: stime ! 初期位置での時刻 [s]
integer, intent(in) :: step ! 計算を行うステップ数
real, intent(in) :: ini_x(:) ! 初期位置 x [m] 複数個設定可
real, intent(in) :: ini_y(size(ini_x)) ! 初期位置 y [m] 複数個設定可
real, intent(in) :: ini_z(size(ini_x)) ! 初期位置 z [m] 複数個設定可
real, intent(in) :: t(:) ! 速度場データのある時刻 [s]
real, intent(in) :: x(:) ! x 方向の座標値 [m]
real, intent(in) :: y(:) ! y 方向の座標値 [m]
real, intent(in) :: z(:) ! z 方向の座標値 [m]
real, intent(in) :: u(size(x),size(y),size(z),size(t)) ! x 方向の風速成分 [m/s]
real, intent(in) :: v(size(x),size(y),size(z),size(t)) ! y 方向の風速成分 [m/s]
real, intent(in) :: w(size(x),size(y),size(z),size(t)) ! z 方向の風速成分 [m/s]
real, intent(inout) :: trajx(step,(size(ini_x))) ! トラジェクトリの位置座標 x [m]
real, intent(inout) :: trajy(step,(size(ini_x))) ! トラジェクトリの位置座標 y [m]
real, intent(inout) :: trajz(step,(size(ini_x))) ! トラジェクトリの位置座標 z [m]
logical, intent(inout) :: FTF(size(ini_x)) ! パーセルが領域外に出たか
! 出たら true,
! 出なければ false.
character(*), intent(in), optional :: opt ! 時間積分スキーム
real, intent(in), optional :: undef ! 未定義
logical, intent(in), optional :: stdopt ! 計算領域外警告除外フラグ
! デフォルトは false (表示する)
character(3) :: option
integer :: nx, ny, nz, nt, nc, i, j, k, l, i_trajt
integer :: defun_i
real :: inter_u(size(x),size(y),size(z))
real :: inter_v(size(x),size(y),size(z))
real :: inter_w(size(x),size(y),size(z))
real :: trajt(step)
real :: defun
logical :: stderr
nx=size(x)
ny=size(y)
nz=size(z)
nt=size(t)
nc=size(ini_x)
! FTF=.false.
trajt=(/((stime+dt*(i-1)),i=1,step)/)
do i=1,nc
trajx(1,i)=ini_x(i)
trajy(1,i)=ini_y(i)
trajz(1,i)=ini_z(i)
end do
!-- 速度場データが stime 以前にない場合, 内挿できないので, エラー
if(trajt(1)<t(1))then
write(*,*) "#### ERROR ####"
write(*,*) "stime must be more than t(1)."
write(*,*) "STOP!!"
stop
end if
if(trajt(step)>t(nt))then
write(*,*) "#### ERROR ####"
write(*,*) "trajt(step) must be less than t(nt)."
write(*,*) "STOP!!"
stop
end if
if(present(opt))then
option=opt
else
option='EU1'
end if
if(present(undef))then
defun=undef
defun_i=int(undef)
else
defun_i=-2147483648
defun=real(defun_i)
end if
do j=1,nc
do i=2,step
trajx(i,j)=defun
trajy(i,j)=defun
trajz(i,j)=defun
end do
end do
if(present(stdopt))then
stderr=stdopt
else
stderr=.false.
end if
do i=1,step-1
!-- 速度場の内挿
if(t(i)==trajt(i).or.t(nt)==trajt(i))then
call interpo_search_1d( t, trajt(i), i_trajt )
if(present(undef))then
do l=1,nz
do k=1,ny
do j=1,nx
if(u(j,k,l,i_trajt)/=undef.and.v(j,k,l,i_trajt)/=undef .and.w(j,k,l,i_trajt)/=undef)then
inter_u(j,k,l)=u(j,k,l,i_trajt)
inter_v(j,k,l)=v(j,k,l,i_trajt)
inter_w(j,k,l)=w(j,k,l,i_trajt)
else
inter_u(j,k,l)=undef
inter_v(j,k,l)=undef
inter_w(j,k,l)=undef
end if
end do
end do
end do
else
do l=1,nz
do k=1,ny
do j=1,nx
inter_u(j,k,l)=u(j,k,l,i_trajt)
inter_v(j,k,l)=v(j,k,l,i_trajt)
inter_w(j,k,l)=w(j,k,l,i_trajt)
end do
end do
end do
end if
else
if(present(undef))then
!$omp parallel default(shared)
!$omp do schedule(dynamic) private(j,k,l)
do l=1,nz
do k=1,ny
do j=1,nx
if(u(j,k,l,i_trajt)/=undef.and.u(j,k,l,i_trajt+1)/=undef.and. v(j,k,l,i_trajt)/=undef.and.v(j,k,l,i_trajt+1)/=undef.and. w(j,k,l,i_trajt)/=undef.and.w(j,k,l,i_trajt+1)/=undef)then
call interpolation_1d( (/t(i_trajt), t(i_trajt+1)/), (/u(j,k,l,i_trajt), u(j,k,l,i_trajt+1)/), trajt(i), inter_u(j,k,l) )
call interpolation_1d( (/t(i_trajt), t(i_trajt+1)/), (/v(j,k,l,i_trajt), v(j,k,l,i_trajt+1)/), trajt(i), inter_v(j,k,l) )
call interpolation_1d( (/t(i_trajt), t(i_trajt+1)/), (/w(j,k,l,i_trajt), w(j,k,l,i_trajt+1)/), trajt(i), inter_w(j,k,l) )
else
inter_u(j,k,l)=undef
inter_v(j,k,l)=undef
inter_w(j,k,l)=undef
end if
end do
end do
end do
!$omp end do
!$omp end parallel
else
!$omp parallel default(shared)
!$omp do schedule(runtime) private(j,k,l)
do l=1,nz
do k=1,ny
do j=1,nx
call interpolation_1d( (/t(i_trajt), t(i_trajt+1)/), (/u(j,k,l,i_trajt), u(j,k,l,i_trajt+1)/), trajt(i), inter_u(j,k,l) )
call interpolation_1d( (/t(i_trajt), t(i_trajt+1)/), (/v(j,k,l,i_trajt), v(j,k,l,i_trajt+1)/), trajt(i), inter_v(j,k,l) )
call interpolation_1d( (/t(i_trajt), t(i_trajt+1)/), (/w(j,k,l,i_trajt), w(j,k,l,i_trajt+1)/), trajt(i), inter_w(j,k,l) )
end do
end do
end do
!$omp end do
!$omp end parallel
end if
end if
if(present(undef))then
!$omp parallel default(shared)
!$omp do schedule(dynamic) private(k)
do k=1,nc
if(FTF(k).eqv..false.)then
call Stream_Line_3d( dt, 2, trajx(i,k), trajy(i,k), trajz(i,k), x, y, z, inter_u, inter_v, inter_w, trajx(i:i+1,k), trajy(i:i+1,k), trajz(i:i+1,k), FTF(k), opt=option, undef=undef, stdopt=stderr )
end if
end do
!$omp end do
!$omp end parallel
else
!$omp parallel default(shared)
!$omp do schedule(dynamic) private(k)
do k=1,nc
if(FTF(k).eqv..false.)then
call Stream_Line_3d( dt, 2, trajx(i,k), trajy(i,k), trajz(i,k), x, y, z, inter_u, inter_v, inter_w, trajx(i:i+1,k), trajy(i:i+1,k), trajz(i:i+1,k), FTF(k), opt=option, stdopt=stderr )
end if
end do
!$omp end do
!$omp end parallel
end if
end do
end subroutine Forward_Traject_3d
| Subroutine : | |||
| dt : | real, intent(in)
| ||
| step : | integer, intent(in)
| ||
| ini_x : | real, intent(in)
| ||
| ini_y : | real, intent(in)
| ||
| x(:) : | real, intent(in)
| ||
| y(:) : | real, intent(in)
| ||
| u(size(x),size(y)) : | real, intent(in)
| ||
| v(size(x),size(y)) : | real, intent(in)
| ||
| trajx(step) : | real, intent(inout)
| ||
| trajy(step) : | real, intent(inout)
| ||
| FTF : | logical, intent(inout)
| ||
| opt : | character(3), intent(in), optional
| ||
| undef : | real, intent(in), optional
| ||
| stdopt : | logical, intent(in), optional
|
流線計算ルーチン 時間積分スキームについては, デフォルトではオイラー法で計算. ‘HO1’= 修正オイラー法(ホインスキーム) ‘ME1’= 改良オイラー法 ‘RK3’= 3 次ルンゲクッタ ‘RK4’= 4 次ルンゲクッタ
subroutine Stream_Line_2d( dt, step, ini_x, ini_y, x, y, u, v, trajx, trajy, FTF, opt, undef, stdopt )
! 流線計算ルーチン
! 時間積分スキームについては, デフォルトではオイラー法で計算.
! '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]
logical, intent(inout) :: FTF ! 各パーセルが領域外に出たか.
! 出て入いれば true,
! 出ていなければ false
character(3), intent(in), optional :: opt ! 時間積分スキーム
real, intent(in), optional :: undef ! 未定義
logical, intent(in), optional :: stdopt ! 計算領域外警告除外フラグ
! デフォルトは false (表示する)
character(3) :: option
!-- 以下, 作業用配列の定義
integer :: i, m, n, id
integer :: nx, ny
integer :: defun_i
real :: k1, k2, k3, k4, l1, l2, l3, l4
real :: x1, y1
real :: inter_p(2) ! 内挿点での位置座標(毎回更新)
real :: inter_v(2) ! 内挿点での速度の値(1=u, 2=v)
real :: defun
logical :: stderr
!-- 汎用性を持たせるために, x, y がどのような範囲をとっても計算可能な仕様にする
!-- そのため, x, y の範囲の要素を限定する
! call (1)
! call (nx)
! call (1)
! call (ny)
FTF=.false.
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
defun_i=int(undef)
else
defun_i=-2147483648
defun=real(defun_i)
end if
if(present(stdopt))then
stderr=stdopt
else
stderr=.false.
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, defun_i )
!-- 速度場が undef であれば, その時点でそれ以降のデータには undef を
!-- 代入し, undef 影響範囲内を通過した旨を伝える.
!-- また, interpo_search において検索点が領域外の場合もそれ以降計算不可能
!-- である旨を通知し, undef があれば undef をなければゼロを代入して抜ける.
if(m==defun_i.or.n==defun_i.or.m==nx.or.n==ny)then
if(stderr.eqv..true.)then
write(*,*) "*********** WARNING ********"
write(*,*) "There is no parcel point in the calculating region."
write(*,*) "Exit!!"
end if
call undef_interpo_1d( step-i+1, trajx(i:step), defun )
call undef_interpo_1d( step-i+1, trajy(i:step), defun )
FTF = .true.
exit
end if
if(present(undef))then
if(u(m,n)==undef.or.u(m,n+1)==undef.or. u(m+1,n)==undef.or.u(m+1,n+1)==undef.or. v(m,n)==undef.or.v(m,n+1)==undef.or. v(m+1,n)==undef.or.v(m+1,n+1)==undef)then
if(stderr.eqv..true.)then
write(*,*) "*********** WARNING ********"
write(*,*) "parcel passed the undef point."
write(*,*) "Exit!!"
end if
call undef_interpo_1d( step-i+1, trajx(i:step), defun )
call undef_interpo_1d( step-i+1, trajy(i:step), defun )
FTF = .true.
exit
end if
end if
!-- 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
if(stderr.eqv..true.)then
write(*,*) "****** WARNING ******"
write(*,*) "This point does not exist in the region."
write(*,*) "Exit.!!"
end if
!-- 計算領域外に出た場合は, 以降のステップにすべて defun を代入
call undef_interpo_1d( step-i+1, trajx(i:step), defun )
call undef_interpo_1d( step-i+1, trajy(i:step), defun )
FTF = .true.
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)/)
do id=1,2
call interpo_search_2d( x, y, inter_p(1), inter_p(2), m, n )
!-- 速度場が undef であれば, その時点でそれ以降のデータには undef を
!-- 代入し, undef 影響範囲内を通過した旨を伝える.
!-- また, interpo_search において検索点が領域外の場合もそれ以降計算不可能
!-- である旨を通知し, undef があれば undef をなければゼロを代入して抜ける.
if(m<=0.or.n<=0.or.m==nx.or.n==ny)then
if(stderr.eqv..true.)then
write(*,*) "*********** WARNING ********"
write(*,*) "There is no parcel point in the calculating region."
write(*,*) "Exit!!"
end if
call undef_interpo_1d( step-i+1, trajx(i:step), defun )
call undef_interpo_1d( step-i+1, trajy(i:step), defun )
FTF = .true.
exit
end if
if(present(undef))then
if(u(m,n)==undef.or.u(m,n+1)==undef.or. u(m+1,n)==undef.or.u(m+1,n+1)==undef.or. v(m,n)==undef.or.v(m,n+1)==undef.or. v(m+1,n)==undef.or.v(m+1,n+1)==undef)then
if(stderr.eqv..true.)then
write(*,*) "*********** WARNING ********"
write(*,*) "parcel passed the undef point."
write(*,*) "Exit!!"
end if
call undef_interpo_1d( step-i+1, trajx(i:step), defun )
call undef_interpo_1d( step-i+1, trajy(i:step), defun )
FTF = .true.
exit
end if
end if
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) )
select case (id)
case (1)
k1=dt*inter_v(1)
l1=dt*inter_v(2)
!-- 一時的な流跡線の位置を計算
x1=trajx(i)+k1
y1=trajy(i)+l1
!-- この点での速度場を内挿
inter_p=(/x1, y1/)
case (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)
end select
end do
!-- 計算した traj が領域内に存在しているか確認
if(x1<x(1).or.x1>x(nx).or.y1<y(1).or.y1>y(ny))then
if(stderr.eqv..true.)then
write(*,*) "****** WARNING ******"
write(*,*) "This point does not exist in the region."
write(*,*) "Exit.!!"
end if
!-- 計算領域外に出た場合は, 以降のステップにすべて defun を代入
call undef_interpo_1d( step-i+1, trajx(i:step), defun )
call undef_interpo_1d( step-i+1, trajy(i:step), defun )
FTF = .true.
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)/)
do id=1,2
call interpo_search_2d( x, y, inter_p(1), inter_p(2), m, n )
!-- 速度場が undef であれば, その時点でそれ以降のデータには undef を
!-- 代入し, undef 影響範囲内を通過した旨を伝える.
!-- また, interpo_search において検索点が領域外の場合もそれ以降計算不可能
!-- である旨を通知し, undef があれば undef をなければゼロを代入して抜ける.
if(m<=0.or.n<=0.or.m==nx.or.n==ny)then
if(stderr.eqv..true.)then
write(*,*) "*********** WARNING ********"
write(*,*) "There is no parcel point in the calculating region."
write(*,*) "Exit!!"
end if
call undef_interpo_1d( step-i+1, trajx(i:step), defun )
call undef_interpo_1d( step-i+1, trajy(i:step), defun )
FTF = .true.
exit
end if
if(present(undef))then
if(u(m,n)==undef.or.u(m,n+1)==undef.or. u(m+1,n)==undef.or.u(m+1,n+1)==undef.or. v(m,n)==undef.or.v(m,n+1)==undef.or. v(m+1,n)==undef.or.v(m+1,n+1)==undef)then
if(stderr.eqv..true.)then
write(*,*) "*********** WARNING ********"
write(*,*) "parcel passed the undef point."
write(*,*) "Exit!!"
end if
call undef_interpo_1d( step-i+1, trajx(i:step), defun )
call undef_interpo_1d( step-i+1, trajy(i:step), defun )
FTF = .true.
exit
end if
end if
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) )
select case (id)
case (1)
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/)
case (2)
k2=dt*inter_v(1)
l2=dt*inter_v(2)
x1=trajx(i)+k2
y1=trajy(i)+l2
end select
end do
!-- 計算した traj が領域内に存在しているか確認
if(x1<x(1).or.x1>x(nx).or.y1<y(1).or.y1>y(ny))then
if(stderr.eqv..true.)then
write(*,*) "****** WARNING ******"
write(*,*) "This point does not exist in the region."
write(*,*) "Exit.!!"
end if
!-- 計算領域外に出た場合は, 以降のステップにすべて defun を代入
call undef_interpo_1d( step-i+1, trajx(i:step), defun )
call undef_interpo_1d( step-i+1, trajy(i:step), defun )
FTF = .true.
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)/)
do id=1,3
call interpo_search_2d( x, y, inter_p(1), inter_p(2), m, n )
!-- 速度場が undef であれば, その時点でそれ以降のデータには undef を
!-- 代入し, undef 影響範囲内を通過した旨を伝える.
!-- また, interpo_search において検索点が領域外の場合もそれ以降計算不可能
!-- である旨を通知し, undef があれば undef をなければゼロを代入して抜ける.
if(m<=0.or.n<=0.or.m==nx.or.n==ny)then
if(stderr.eqv..true.)then
write(*,*) "*********** WARNING ********"
write(*,*) "There is no parcel point in the calculating region."
write(*,*) "Exit!!"
end if
call undef_interpo_1d( step-i+1, trajx(i:step), defun )
call undef_interpo_1d( step-i+1, trajy(i:step), defun )
FTF = .true.
exit
end if
if(present(undef))then
if(u(m,n)==undef.or.u(m,n+1)==undef.or. u(m+1,n)==undef.or.u(m+1,n+1)==undef.or. v(m,n)==undef.or.v(m,n+1)==undef.or. v(m+1,n)==undef.or.v(m+1,n+1)==undef)then
if(stderr.eqv..true.)then
write(*,*) "*********** WARNING ********"
write(*,*) "parcel passed the undef point."
write(*,*) "Exit!!"
end if
call undef_interpo_1d( step-i+1, trajx(i:step), defun )
call undef_interpo_1d( step-i+1, trajy(i:step), defun )
FTF = .true.
exit
end if
end if
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) )
select case (id)
case (1)
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/)
case (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/)
case (3)
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)
end select
end do
!-- 計算した traj が領域内に存在しているか確認
if(x1<x(1).or.x1>x(nx).or.y1<y(1).or.y1>y(ny))then
if(stderr.eqv..true.)then
write(*,*) "****** WARNING ******"
write(*,*) "This point does not exist in the region."
write(*,*) "Exit.!!"
end if
!-- 計算領域外に出た場合は, 以降のステップにすべて defun を代入
call undef_interpo_1d( step-i+1, trajx(i:step), defun )
call undef_interpo_1d( step-i+1, trajy(i:step), defun )
FTF = .true.
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)/)
do id=1,4
call interpo_search_2d( x, y, inter_p(1), inter_p(2), m, n )
!-- 速度場が undef であれば, その時点でそれ以降のデータには undef を
!-- 代入し, undef 影響範囲内を通過した旨を伝える.
!-- また, interpo_search において検索点が領域外の場合もそれ以降計算不可能
!-- である旨を通知し, undef があれば undef をなければゼロを代入して抜ける.
if(m<=0.or.n<=0.or.m==nx.or.n==ny)then
if(stderr.eqv..true.)then
write(*,*) "*********** WARNING ********"
write(*,*) "There is no parcel point in the calculating region."
write(*,*) "Exit!!"
end if
call undef_interpo_1d( step-i+1, trajx(i:step), defun )
call undef_interpo_1d( step-i+1, trajy(i:step), defun )
FTF = .true.
exit
end if
if(present(undef))then
if(u(m,n)==undef.or.u(m,n+1)==undef.or. u(m+1,n)==undef.or.u(m+1,n+1)==undef.or. v(m,n)==undef.or.v(m,n+1)==undef.or. v(m+1,n)==undef.or.v(m+1,n+1)==undef)then
if(stderr.eqv..true.)then
write(*,*) "*********** WARNING ********"
write(*,*) "parcel passed the undef point."
write(*,*) "Exit!!"
end if
call undef_interpo_1d( step-i+1, trajx(i:step), defun )
call undef_interpo_1d( step-i+1, trajy(i:step), defun )
FTF = .true.
exit
end if
end if
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) )
select case (id)
case (1)
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/)
case (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/)
case (3)
k3=dt*inter_v(1)
l3=dt*inter_v(2)
x1=trajx(i)+k3
y1=trajy(i)+l3
!-- この点での速度場を内挿
inter_p=(/x1, y1/)
case (4)
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)
end select
end do
!-- 計算した traj が領域内に存在しているか確認
if(x1<x(1).or.x1>x(nx).or.y1<y(1).or.y1>y(ny))then
if(stderr.eqv..true.)then
write(*,*) "****** WARNING ******"
write(*,*) "This point does not exist in the region."
write(*,*) "Exit.!!"
end if
!-- 計算領域外に出た場合は, 以降のステップにすべて defun を代入
call undef_interpo_1d( step-i+1, trajx(i:step), defun )
call undef_interpo_1d( step-i+1, trajy(i:step), defun )
FTF = .true.
exit
end if
trajx(i+1)=x1
trajy(i+1)=y1
end do
end select
end subroutine Stream_Line_2d
| Subroutine : | |||
| dt : | real, intent(in)
| ||
| step : | integer, intent(in)
| ||
| ini_x : | real, intent(in)
| ||
| ini_y : | real, intent(in)
| ||
| ini_z : | real, intent(in)
| ||
| x(:) : | real, intent(in)
| ||
| y(:) : | real, intent(in)
| ||
| z(:) : | real, intent(in)
| ||
| u(size(x),size(y),size(z)) : | real, intent(in)
| ||
| v(size(x),size(y),size(z)) : | real, intent(in)
| ||
| w(size(x),size(y),size(z)) : | real, intent(in)
| ||
| trajx(step) : | real, intent(inout)
| ||
| trajy(step) : | real, intent(inout)
| ||
| trajz(step) : | real, intent(inout)
| ||
| FTF : | logical, intent(inout)
| ||
| opt : | character(*), intent(in), optional
| ||
| undef : | real, intent(in), optional
| ||
| stdopt : | logical, intent(in), optional
|
流線計算ルーチン 時間積分スキームについては, デフォルトではオイラー法で計算. ‘HO1’= 修正オイラー法(ホインスキーム) ‘ME1’= 改良オイラー法 ‘RK3’= 3 次ルンゲクッタ ‘RK4’= 4 次ルンゲクッタ
subroutine Stream_Line_3d( dt, step, ini_x, ini_y, ini_z, x, y, z, u, v, w, trajx, trajy, trajz, FTF, opt, undef, stdopt )
! 流線計算ルーチン
! 時間積分スキームについては, デフォルトではオイラー法で計算.
! '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]
logical, intent(inout) :: FTF ! 各パーセルが領域外に出たか.
! 出て入いれば true,
! 出ていなければ false
character(*), intent(in), optional :: opt ! 時間積分スキーム
real, intent(in), optional :: undef ! 未定義
logical, intent(in), optional :: stdopt ! 計算領域外警告除外フラグ
! デフォルトは false (表示する)
character(3) :: option
!-- 以下, 作業用配列の定義
integer :: i, l, m, n, nx, ny, nz, id
integer :: defun_i
real :: k1, k2, k3, k4, l1, l2, l3, l4
real :: x1, y1, z1, n1, n2, n3, n4
real :: inter_p(3) ! 内挿点での位置座標(毎回更新)
real :: inter_v(3) ! 内挿点での速度の値(1=u, 2=v)
real :: defun
logical :: stderr
!-- 汎用性を持たせるために, x, y がどのような範囲をとっても計算可能な仕様にする
!-- そのため, x, y の範囲の要素を限定する
! call (1)
! call (nx)
! call (1)
! call (ny)
FTF=.false.
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
if(present(undef))then
defun=undef
defun_i=int(undef)
else
defun_i=-2147483648
defun=real(defun_i)
end if
if(present(stdopt))then
stderr=stdopt
else
stderr=.false.
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 )
!-- 速度場が undef であれば, その時点でそれ以降のデータには undef を
!-- 代入し, undef 影響範囲内を通過した旨を伝える.
!-- また, interpo_search において検索点が領域外の場合もそれ以降計算不可能
!-- である旨を通知し, undef があれば undef をなければゼロを代入して抜ける.
if(m<=0.or.n<=0.or.l<=0.or.m==nx.or.n==ny.or.l==nz)then
if(stderr.eqv..true.)then
write(*,*) "*********** WARNING ********"
write(*,*) "There is no parcel point in the calculating region."
write(*,*) "Exit!!"
end if
call undef_interpo_1d( step-i+1, trajx(i:step), defun )
call undef_interpo_1d( step-i+1, trajy(i:step), defun )
call undef_interpo_1d( step-i+1, trajz(i:step), defun )
FTF = .true.
exit
end if
if(present(undef))then
if(u(m+1,n+1,l+1)==undef.or.u(m,n,l)==undef.or. u(m,n+1,l)==undef.or.u(m,n,l+1)==undef.or. u(m+1,n,l)==undef.or.u(m+1,n,l+1)==undef.or. u(m+1,n+1,l)==undef.or.u(m,n+1,l+1)==undef.or. v(m+1,n+1,l+1)==undef.or.v(m,n,l)==undef.or. v(m,n+1,l)==undef.or.v(m,n,l+1)==undef.or. v(m+1,n,l)==undef.or.v(m+1,n,l+1)==undef.or. v(m+1,n+1,l)==undef.or.v(m,n+1,l+1)==undef.or. w(m+1,n+1,l+1)==undef.or.w(m,n,l)==undef.or. w(m,n+1,l)==undef.or.w(m,n,l+1)==undef.or. w(m+1,n,l)==undef.or.w(m+1,n,l+1)==undef.or. w(m+1,n+1,l)==undef.or.w(m,n+1,l+1)==undef)then
if(stderr.eqv..true.)then
write(*,*) "*********** WARNING ********"
write(*,*) "parcel passed the undef point."
write(*,*) "Exit!!"
end if
call undef_interpo_1d( step-i+1, trajx(i:step), defun )
call undef_interpo_1d( step-i+1, trajy(i:step), defun )
call undef_interpo_1d( step-i+1, trajz(i:step), defun )
FTF = .true.
exit
end if
end if
!-- 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
if(stderr.eqv..true.)then
write(*,*) "****** WARNING ******"
write(*,*) "This point does not exist in the region."
write(*,*) "Exit.!!"
end if
!-- 計算領域外に出た場合は, 以降のステップにすべて defun を代入
call undef_interpo_1d( step-i+1, trajx(i:step), defun )
call undef_interpo_1d( step-i+1, trajy(i:step), defun )
call undef_interpo_1d( step-i+1, trajz(i:step), defun )
FTF = .true.
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)/)
do id=1,2
call interpo_search_3d( x, y, z, inter_p(1), inter_p(2), inter_p(3), m, n, l )
!-- 速度場が undef であれば, その時点でそれ以降のデータには undef を
!-- 代入し, undef 影響範囲内を通過した旨を伝える.
!-- また, interpo_search において検索点が領域外の場合もそれ以降計算不可能
!-- である旨を通知し, undef があれば undef をなければゼロを代入して抜ける.
if(m<=0.or.n<=0.or.l<=0.or.m==nx.or.n==ny.or.l==nz)then
if(stderr.eqv..true.)then
write(*,*) "*********** WARNING ********"
write(*,*) "There is no parcel point in the calculating region."
write(*,*) "Exit!!"
end if
call undef_interpo_1d( step-i+1, trajx(i:step), defun )
call undef_interpo_1d( step-i+1, trajy(i:step), defun )
call undef_interpo_1d( step-i+1, trajz(i:step), defun )
FTF = .true.
exit
end if
if(present(undef))then
if(u(m+1,n+1,l+1)==undef.or.u(m,n,l)==undef.or. u(m,n+1,l)==undef.or.u(m,n,l+1)==undef.or. u(m+1,n,l)==undef.or.u(m+1,n,l+1)==undef.or. u(m+1,n+1,l)==undef.or.u(m,n+1,l+1)==undef.or. v(m+1,n+1,l+1)==undef.or.v(m,n,l)==undef.or. v(m,n+1,l)==undef.or.v(m,n,l+1)==undef.or. v(m+1,n,l)==undef.or.v(m+1,n,l+1)==undef.or. v(m+1,n+1,l)==undef.or.v(m,n+1,l+1)==undef.or. w(m+1,n+1,l+1)==undef.or.w(m,n,l)==undef.or. w(m,n+1,l)==undef.or.w(m,n,l+1)==undef.or. w(m+1,n,l)==undef.or.w(m+1,n,l+1)==undef.or. w(m+1,n+1,l)==undef.or.w(m,n+1,l+1)==undef)then
if(stderr.eqv..true.)then
write(*,*) "*********** WARNING ********"
write(*,*) "parcel passed the undef point."
write(*,*) "Exit!!"
end if
call undef_interpo_1d( step-i+1, trajx(i:step), defun )
call undef_interpo_1d( step-i+1, trajy(i:step), defun )
call undef_interpo_1d( step-i+1, trajz(i:step), defun )
FTF = .true.
exit
end if
end if
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) )
select case (id)
case (1)
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/)
case (2)
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)
end select
end do
!-- 計算した 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
if(stderr.eqv..true.)then
write(*,*) "****** WARNING ******"
write(*,*) "This point does not exist in the region."
write(*,*) "Exit.!!"
end if
!-- 計算領域外に出た場合は, 以降のステップにすべて defun を代入
call undef_interpo_1d( step-i+1, trajx(i:step), defun )
call undef_interpo_1d( step-i+1, trajy(i:step), defun )
call undef_interpo_1d( step-i+1, trajz(i:step), defun )
FTF = .true.
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)/)
do id=1,2
call interpo_search_3d( x, y, z, inter_p(1), inter_p(2), inter_p(3), m, n, l )
!-- 速度場が undef であれば, その時点でそれ以降のデータには undef を
!-- 代入し, undef 影響範囲内を通過した旨を伝える.
!-- また, interpo_search において検索点が領域外の場合もそれ以降計算不可能
!-- である旨を通知し, undef があれば undef をなければゼロを代入して抜ける.
if(m<=0.or.n<=0.or.l<=0.or.m==nx.or.n==ny.or.l==nz)then
if(stderr.eqv..true.)then
write(*,*) "*********** WARNING ********"
write(*,*) "There is no parcel point in the calculating region."
write(*,*) "Exit!!"
end if
call undef_interpo_1d( step-i+1, trajx(i:step), defun )
call undef_interpo_1d( step-i+1, trajy(i:step), defun )
call undef_interpo_1d( step-i+1, trajz(i:step), defun )
FTF = .true.
exit
end if
if(present(undef))then
if(u(m+1,n+1,l+1)==undef.or.u(m,n,l)==undef.or. u(m,n+1,l)==undef.or.u(m,n,l+1)==undef.or. u(m+1,n,l)==undef.or.u(m+1,n,l+1)==undef.or. u(m+1,n+1,l)==undef.or.u(m,n+1,l+1)==undef.or. v(m+1,n+1,l+1)==undef.or.v(m,n,l)==undef.or. v(m,n+1,l)==undef.or.v(m,n,l+1)==undef.or. v(m+1,n,l)==undef.or.v(m+1,n,l+1)==undef.or. v(m+1,n+1,l)==undef.or.v(m,n+1,l+1)==undef.or. w(m+1,n+1,l+1)==undef.or.w(m,n,l)==undef.or. w(m,n+1,l)==undef.or.w(m,n,l+1)==undef.or. w(m+1,n,l)==undef.or.w(m+1,n,l+1)==undef.or. w(m+1,n+1,l)==undef.or.w(m,n+1,l+1)==undef)then
if(stderr.eqv..true.)then
write(*,*) "*********** WARNING ********"
write(*,*) "parcel passed the undef point."
write(*,*) "Exit!!"
end if
call undef_interpo_1d( step-i+1, trajx(i:step), defun )
call undef_interpo_1d( step-i+1, trajy(i:step), defun )
call undef_interpo_1d( step-i+1, trajz(i:step), defun )
FTF = .true.
exit
end if
end if
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) )
select case (id)
case (1)
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/)
case (2)
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
end select
end do
!-- 計算した 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
if(stderr.eqv..true.)then
write(*,*) "****** WARNING ******"
write(*,*) "This point does not exist in the region."
write(*,*) "Exit.!!"
end if
!-- 計算領域外に出た場合は, 以降のステップにすべて defun を代入
call undef_interpo_1d( step-i+1, trajx(i:step), defun )
call undef_interpo_1d( step-i+1, trajy(i:step), defun )
call undef_interpo_1d( step-i+1, trajz(i:step), defun )
FTF = .true.
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)/)
do id=1,3
call interpo_search_3d( x, y, z, inter_p(1), inter_p(2), inter_p(3), m, n, l )
!-- 速度場が undef であれば, その時点でそれ以降のデータには undef を
!-- 代入し, undef 影響範囲内を通過した旨を伝える.
!-- また, interpo_search において検索点が領域外の場合もそれ以降計算不可能
!-- である旨を通知し, undef があれば undef をなければゼロを代入して抜ける.
if(m<=0.or.n<=0.or.l<=0.or.m==nx.or.n==ny.or.l==nz)then
if(stderr.eqv..true.)then
write(*,*) "*********** WARNING ********"
write(*,*) "There is no parcel point in the calculating region."
write(*,*) "Exit!!"
end if
call undef_interpo_1d( step-i+1, trajx(i:step), defun )
call undef_interpo_1d( step-i+1, trajy(i:step), defun )
call undef_interpo_1d( step-i+1, trajz(i:step), defun )
FTF = .true.
exit
end if
if(present(undef))then
if(u(m+1,n+1,l+1)==undef.or.u(m,n,l)==undef.or. u(m,n+1,l)==undef.or.u(m,n,l+1)==undef.or. u(m+1,n,l)==undef.or.u(m+1,n,l+1)==undef.or. u(m+1,n+1,l)==undef.or.u(m,n+1,l+1)==undef.or. v(m+1,n+1,l+1)==undef.or.v(m,n,l)==undef.or. v(m,n+1,l)==undef.or.v(m,n,l+1)==undef.or. v(m+1,n,l)==undef.or.v(m+1,n,l+1)==undef.or. v(m+1,n+1,l)==undef.or.v(m,n+1,l+1)==undef.or. w(m+1,n+1,l+1)==undef.or.w(m,n,l)==undef.or. w(m,n+1,l)==undef.or.w(m,n,l+1)==undef.or. w(m+1,n,l)==undef.or.w(m+1,n,l+1)==undef.or. w(m+1,n+1,l)==undef.or.w(m,n+1,l+1)==undef)then
if(stderr.eqv..true.)then
write(*,*) "*********** WARNING ********"
write(*,*) "parcel passed the undef point."
write(*,*) "Exit!!"
end if
call undef_interpo_1d( step-i+1, trajx(i:step), defun )
call undef_interpo_1d( step-i+1, trajy(i:step), defun )
call undef_interpo_1d( step-i+1, trajz(i:step), defun )
FTF = .true.
exit
end if
end if
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) )
select case (id)
case (1)
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/)
case (2)
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/)
case (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)
end select
end do
!-- 計算した 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
if(stderr.eqv..true.)then
write(*,*) "****** WARNING ******"
write(*,*) "This point does not exist in the region."
write(*,*) "Exit.!!"
end if
!-- 計算領域外に出た場合は, 以降のステップにすべて defun を代入
call undef_interpo_1d( step-i+1, trajx(i:step), defun )
call undef_interpo_1d( step-i+1, trajy(i:step), defun )
call undef_interpo_1d( step-i+1, trajz(i:step), defun )
FTF = .true.
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)/)
do id=1,4
call interpo_search_3d( x, y, z, inter_p(1), inter_p(2), inter_p(3), m, n, l )
!-- 速度場が undef であれば, その時点でそれ以降のデータには undef を
!-- 代入し, undef 影響範囲内を通過した旨を伝える.
!-- また, interpo_search において検索点が領域外の場合もそれ以降計算不可能
!-- である旨を通知し, undef があれば undef をなければゼロを代入して抜ける.
if(m<=0.or.n<=0.or.l<=0.or.m==nx.or.n==ny.or.l==nz)then
if(stderr.eqv..true.)then
write(*,*) "*********** WARNING ********"
write(*,*) "There is no parcel point in the calculating region."
write(*,*) "Exit!!"
end if
call undef_interpo_1d( step-i+1, trajx(i:step), defun )
call undef_interpo_1d( step-i+1, trajy(i:step), defun )
call undef_interpo_1d( step-i+1, trajz(i:step), defun )
FTF = .true.
exit
end if
if(present(undef))then
if(u(m+1,n+1,l+1)==undef.or.u(m,n,l)==undef.or. u(m,n+1,l)==undef.or.u(m,n,l+1)==undef.or. u(m+1,n,l)==undef.or.u(m+1,n,l+1)==undef.or. u(m+1,n+1,l)==undef.or.u(m,n+1,l+1)==undef.or. v(m+1,n+1,l+1)==undef.or.v(m,n,l)==undef.or. v(m,n+1,l)==undef.or.v(m,n,l+1)==undef.or. v(m+1,n,l)==undef.or.v(m+1,n,l+1)==undef.or. v(m+1,n+1,l)==undef.or.v(m,n+1,l+1)==undef.or. w(m+1,n+1,l+1)==undef.or.w(m,n,l)==undef.or. w(m,n+1,l)==undef.or.w(m,n,l+1)==undef.or. w(m+1,n,l)==undef.or.w(m+1,n,l+1)==undef.or. w(m+1,n+1,l)==undef.or.w(m,n+1,l+1)==undef)then
if(stderr.eqv..true.)then
write(*,*) "*********** WARNING ********"
write(*,*) "parcel passed the undef point."
write(*,*) "Exit!!"
end if
call undef_interpo_1d( step-i+1, trajx(i:step), defun )
call undef_interpo_1d( step-i+1, trajy(i:step), defun )
call undef_interpo_1d( step-i+1, trajz(i:step), defun )
FTF = .true.
exit
end if
end if
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) )
select case (id)
case (1)
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/)
case (2)
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/)
case (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/)
case (4)
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)
end select
end do
!-- 計算した 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
if(stderr.eqv..true.)then
write(*,*) "****** WARNING ******"
write(*,*) "This point does not exist in the region."
write(*,*) "Exit.!!"
end if
!-- 計算領域外に出た場合は, 以降のステップにすべて defun を代入
call undef_interpo_1d( step-i+1, trajx(i:step), defun )
call undef_interpo_1d( step-i+1, trajy(i:step), defun )
call undef_interpo_1d( step-i+1, trajz(i:step), defun )
FTF = .true.
exit
end if
trajx(i+1)=x1
trajy(i+1)=y1
trajz(i+1)=z1
end do
end select
end subroutine Stream_Line_3d