Class analy
In: analy.f90

計算結果データ用解析モジュール

Methods

Public Instance methods

Subroutine :
xp :real, intent(in)
: 回転軸の x 位置座標
yp :real, intent(in)
: 回転軸の y 位置座標
x(:) :real, intent(in)
: x 方向の位置座標
y(:) :real, intent(in)
: y 方向の位置座標
u(size(x),size(y)) :real, intent(in)
: 位置 i,j での風速の 1 成分
v(size(x),size(y)) :real, intent(in)
: 位置 i,j での風速の 1 成分
f(size(x),size(y)) :real, intent(in)
: 位置 i,j でのコリオリパラメータ
mome(size(x),size(y)) :real, intent(inout)
: 回転軸まわりの相対角運動量

任意の点まわりの絶対角運動量を計算するルーチン

主目的は台風の中心軸を中心に鉛直軸まわりの角運動量を計算することを目的とする.

$$M=rv+dfrac{fr^2}{2} ,quad r=中心軸からの距離, quad v=風速の同位角成分$$

位置と風速に緯度の変換を与えれば, 全球での自転軸まわりの角運動量を 計算することも可能.

ベクトルの外積計算ルーチン vec_prod を用いることで極座標でも計算可能.

[Source]

subroutine Aangular_moment(xp,yp,x,y,u,v,f,mome)
! 任意の点まわりの絶対角運動量を計算するルーチン
!
! 主目的は台風の中心軸を中心に鉛直軸まわりの角運動量を計算することを目的とする.
!
! $$M=rv+\dfrac{fr^2}{2} ,\quad r=中心軸からの距離, \quad v=風速の同位角成分$$
!
! 位置と風速に緯度の変換を与えれば, 全球での自転軸まわりの角運動量を
! 計算することも可能.
!
! ベクトルの外積計算ルーチン vec_prod を用いることで極座標でも計算可能.
  implicit none
  real, intent(in) :: x(:)  ! x 方向の位置座標
  real, intent(in) :: y(:)  ! y 方向の位置座標
  real, intent(in) :: xp  ! 回転軸の x 位置座標
  real, intent(in) :: yp  ! 回転軸の y 位置座標
  real, intent(in) :: u(size(x),size(y))  ! 位置 i,j での風速の 1 成分
  real, intent(in) :: v(size(x),size(y))  ! 位置 i,j での風速の 1 成分
  real, intent(in) :: f(size(x),size(y))  ! 位置 i,j でのコリオリパラメータ
  real, intent(inout) :: mome(size(x),size(y))  ! 回転軸まわりの相対角運動量
  real :: xxp(size(x),size(y),1)  ! x,y,xp,yp から計算した軸中心からの位置ベクトル x 成分
  real :: yyp(size(x),size(y),1)  ! x,y,xp,yp から計算した軸中心からの位置ベクトル y 成分
  integer :: i, j, nx, ny
  real :: tmp(size(x),size(y),1), rp(size(x),size(y)), tmp1(1)

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

!$omp parallel do shared(xxp,yyp,x,y,xp,yp) private(i,j)
  do j=1,ny
     do i=1,nx
        xxp(i,j,1)=x(i)-xp
        yyp(i,j,1)=y(j)-yp
     end do
  end do
!$omp end parallel do

  tmp=0.0
  call vec_prod(xxp,yyp,tmp,u,v,tmp,tmp,tmp,mome)
  call radius(xp,yp,0.0,x,y,tmp1,rp)

!$omp parallel do shared(mome,f,rp) private(i,j)
  do j=1,ny
     do i=1,nx
        mome(i,j)=mome(i,j)+0.5*f(i,j)*rp(i,j)**2
     end do
  end do
!$omp end parallel do

end subroutine Aangular_moment
Subroutine :
xp :real, intent(in)
: 回転軸の x 位置座標
yp :real, intent(in)
: 回転軸の y 位置座標
x(:) :real, intent(in)
: x 方向の位置座標
y(:) :real, intent(in)
: y 方向の位置座標
u(size(x),size(y)) :real, intent(in)
: 位置 i,j での風速の 1 成分
v(size(x),size(y)) :real, intent(in)
: 位置 i,j での風速の 1 成分
mome(size(x),size(y)) :real, intent(inout)
: 回転軸まわりの相対角運動量

任意の点まわりの相対角運動量を計算するルーチン

本当は 3 次元ベクトルで計算されるが, 気象学では 3 次元量はあまり需要がない であろうという判断から, ある回転軸まわりの角運動量成分のみを 計算することにしている.

主目的は台風の中心軸を中心に鉛直軸まわりの角運動量を計算することを目的とする.

$$M=rv,quad r=中心軸からの距離, quad v=風速の同位角成分$$

位置と風速に緯度の変換を与えれば, 全球での自転軸まわりの角運動量を 計算することも可能. ベクトルの外積計算ルーチン vec_prod を用いることで極座標でも計算可能.

[Source]

subroutine Rangular_moment(xp,yp,x,y,u,v,mome)
! 任意の点まわりの相対角運動量を計算するルーチン
!
! 本当は 3 次元ベクトルで計算されるが, 気象学では 3 次元量はあまり需要がない
! であろうという判断から, ある回転軸まわりの角運動量成分のみを
! 計算することにしている.
!
! 主目的は台風の中心軸を中心に鉛直軸まわりの角運動量を計算することを目的とする.
!
! $$M=rv,\quad r=中心軸からの距離, \quad v=風速の同位角成分$$
!
! 位置と風速に緯度の変換を与えれば, 全球での自転軸まわりの角運動量を
! 計算することも可能.
! ベクトルの外積計算ルーチン vec_prod を用いることで極座標でも計算可能.
  implicit none
  real, intent(in) :: x(:)  ! x 方向の位置座標
  real, intent(in) :: y(:)  ! y 方向の位置座標
  real, intent(in) :: xp  ! 回転軸の x 位置座標
  real, intent(in) :: yp  ! 回転軸の y 位置座標
  real, intent(in) :: u(size(x),size(y))  ! 位置 i,j での風速の 1 成分
  real, intent(in) :: v(size(x),size(y))  ! 位置 i,j での風速の 1 成分
  real, intent(inout) :: mome(size(x),size(y))  ! 回転軸まわりの相対角運動量
  real :: xxp(size(x),size(y),1)  ! x,y,xp,yp から計算した軸中心からの位置ベクトル x 成分
  real :: yyp(size(x),size(y),1)  ! x,y,xp,yp から計算した軸中心からの位置ベクトル y 成分
  integer :: i, j, nx, ny
  real :: tmp(size(x),size(y),1)

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

!$omp parallel do shared(xxp,yyp,x,y,xp,yp) private(i,j)
  do j=1,ny
     do i=1,nx
        xxp(i,j,1)=x(i)-xp
        yyp(i,j,1)=y(j)-yp
     end do
  end do
!$omp end parallel do

  tmp=0.0
  call vec_prod(xxp,yyp,tmp,u,v,tmp,tmp,tmp,mome)

end subroutine Rangular_moment
Subroutine :
x(:,:,:) :real, intent(in)
: x 方向のベクトル成分
y(size(x,1),size(x,2),size(x,3)) :real, intent(in)
: y 方向のベクトル成分
z(size(x,1),size(x,2),size(x,3)) :real, intent(in)
: z 方向のベクトル成分
dis(size(x,1),size(x,2),size(x,3)) :real, intent(inout)
: 各点での絶対値ベクトル

3 次元ベクトルの絶対値を計算するルーチン 配列数を調整することにより, 2 次元での計算も可能.

[Source]

subroutine abst(x,y,z,dis)  ! 3 次元ベクトルの絶対値を計算するルーチン
  ! 配列数を調整することにより, 2 次元での計算も可能.
  implicit none
  real, intent(in) :: x(:,:,:)  ! x 方向のベクトル成分
  real, intent(in) :: y(size(x,1),size(x,2),size(x,3))  ! y 方向のベクトル成分
  real, intent(in) :: z(size(x,1),size(x,2),size(x,3))  ! z 方向のベクトル成分
  real, intent(inout) :: dis(size(x,1),size(x,2),size(x,3))  ! 各点での絶対値ベクトル
  integer :: i, j, k, nx, ny, nz

  nx=size(x,1)
  ny=size(x,2)
  nz=size(x,3)

!$omp parallel do shared(dis,x,y,z) private(i,j,k)
  do k=1,nz
     do j=1,ny
        do i=1,nx
           dis(i,j,k)=sqrt(x(i,j,k)**2+y(i,j,k)**2+z(i,j,k)**2)
        end do
     end do
  end do
!$omp end parallel do

end subroutine abst
Subroutine :
x(:) :real, intent(in)
: x 方向の空間座標 [m]
y(:) :real, intent(in)
: y 方向の空間座標 [m]
u(size(x),size(y)) :real, intent(in)
: x に対応する方向の 2 次元ベクトル成分
v(size(x),size(y)) :real, intent(in)
: y に対応する方向の 2 次元ベクトル成分
val(size(x),size(y)) :real, intent(inout)
: 2次元渦度
undeff :real, intent(in), optional
hx(size(x),size(y)) :real, intent(in), optional
: x 方向のスケール因子
hy(size(x),size(y)) :real, intent(in), optional
: y 方向のスケール因子
ord :logical, intent(in), optional
: 微分の順番を入れ替えるオプション. true なら, 入れ替えない, false なら, 入れ替える. デフォルトは true, du/dz-dw/dx を計算するときのみ用いる.

2 次元渦度計算ルーチン

x, y は配列の次元の若い順に必ず並べること.

u, v は偶置換の向きに配置すれば, 任意の2次元渦度が計算可能 ただし, du/dz-dw/dx を計算するときのみ, (x,z,u,w) の順番で, ord オプション false.

$frac{partial v}{partial x} -frac{partial u}{partial y} $ を 2 次の中央差分近似で書き換えると, 点 $(i,j)$ での発散は $frac{v_{i,j+1}-v_{i,j-1}}{2dx} -frac{u_{i+1,j}-u_{i-1,j}}{2dy} $ とできる. これを用いて2次元発散を計算. データ点が足りない隅の領域では, 1 次の差分近似で計算するので, 少し精度が 落ちる. 基本的に等間隔格子のデータのみしか計算できないので, 不等間隔のデータは 等間隔に内挿なりを行ってからこのルーチンを呼び出す必要あり. 将来, 任意の格子座標系に対応するため, 後ろにスケーリングファクターオプションを 追加したが, 現在は使用していない.

[Source]

subroutine curl( x, y, u, v, val, undeff, hx, hy, ord )
! 2 次元渦度計算ルーチン
!
! x, y は配列の次元の若い順に必ず並べること.
!
! u, v は偶置換の向きに配置すれば, 任意の2次元渦度が計算可能
! ただし, du/dz-dw/dx を計算するときのみ, (x,z,u,w) の順番で, ord オプション false.
!
! $\frac{\partial v}{\partial x} -\frac{\partial u}{\partial y} $ を
! 2 次の中央差分近似で書き換えると, 点 $(i,j)$ での発散は
! $\frac{v_{i,j+1}-v_{i,j-1}}{2dx} -\frac{u_{i+1,j}-u_{i-1,j}}{2dy} $
! とできる. これを用いて2次元発散を計算.
! データ点が足りない隅の領域では, 1 次の差分近似で計算するので, 少し精度が
! 落ちる.
! 基本的に等間隔格子のデータのみしか計算できないので, 不等間隔のデータは
! 等間隔に内挿なりを行ってからこのルーチンを呼び出す必要あり.
! 将来, 任意の格子座標系に対応するため, 後ろにスケーリングファクターオプションを
! 追加したが, 現在は使用していない.
  implicit none
  real, intent(in) :: x(:)  ! x 方向の空間座標 [m]
  real, intent(in) :: y(:)  ! y 方向の空間座標 [m]
  real, intent(in) :: u(size(x),size(y))  ! x に対応する方向の 2 次元ベクトル成分
  real, intent(in) :: v(size(x),size(y))  ! y に対応する方向の 2 次元ベクトル成分
  real, intent(inout) :: val(size(x),size(y))  ! 2次元渦度
  real, intent(in), optional :: undeff
  real, intent(in), optional :: hx(size(x),size(y))  ! x 方向のスケール因子
  real, intent(in), optional :: hy(size(x),size(y))  ! y 方向のスケール因子
  logical, intent(in),  optional :: ord  ! 微分の順番を入れ替えるオプション.
                 ! true なら, 入れ替えない, false なら, 入れ替える.
                 ! デフォルトは true, du/dz-dw/dx を計算するときのみ用いる.
  integer :: i   ! イタレーション用添字
  integer :: j   ! イタレーション用添字
  integer :: nx  ! 空間配列要素数 1 次元目
  integer :: ny  ! 空間配列要素数 2 次元目
  real :: dx  ! 1 次元目を微分する格子間隔 [m]
  real :: dy  ! 2 次元目を微分する格子間隔 [m]
  logical :: order

  if(present(ord))then
     order=ord
  else
     order=.true.
  end if

  dx=x(2)-x(1)
  dy=y(2)-y(1)
  nx=size(x)
  ny=size(y)

  if(present(undeff))then
     do j=2,ny-1
        do i=2,nx-1
           if(v(i+1,j)==undeff.or.v(i-1,j)==undeff.or.u(i,j+1)==undeff.or. u(i,j-1)==undeff)then
              val(i,j)=undeff
           else
              val(i,j)=0.5*((v(i+1,j)-v(i-1,j))/dx-(u(i,j+1)-u(i,j-1))/dy)
           end if
        end do
     end do

!-- 領域の端 ---
!-- y 方向の両端 ---
     do i=2,nx-1
        if(v(i+1,1)==undeff.or.v(i-1,1)==undeff.or.u(i,2)==undeff.or.u(i,1)==undeff)then
           val(i,1)=undeff
        else
           val(i,1)=0.5*(v(i+1,1)-v(i-1,1))/dx-(u(i,2)-u(i,1))/dy
        end if

        if(v(i+1,ny)==undeff.or.v(i-1,ny)==undeff.or.u(i,ny)==undeff.or. u(i,ny-1)==undeff)then
           val(i,ny)=undeff
        else
           val(i,ny)=0.5*(v(i+1,ny)-v(i-1,ny))/dx-(u(i,ny)-u(i,ny-1))/dy
        end if
     end do
!-- x 方向の両端 ---
     do j=2,ny-1
        if(v(2,j)==undeff.or.v(1,j)==undeff.or.u(1,j+1)==undeff.or.u(1,j-1)==undeff)then
           val(1,j)=undeff
        else
           val(1,j)=(v(2,j)-v(1,j))/dx-0.5*((u(1,j+1)-u(1,j-1))/dy)
        end if

        if(v(nx,j)==undeff.or.v(nx-1,j)==undeff.or.u(nx,j+1)==undeff.or. u(nx,j-1)==undeff)then
           val(nx,j)=undeff
        else
           val(nx,j)=(v(nx,j)-v(nx-1,j))/dx-0.5*((u(nx,j+1)-u(nx,j-1))/dy)
        end if
     end do
!-- 4 隅 ---
     if(v(2,1)==undeff.or.v(1,1)==undeff.or.u(1,2)==undeff.or.u(1,1)==undeff)then
        val(1,1)=undeff
     else
        val(1,1)=(v(2,1)-v(1,1))/dx-(u(1,2)-u(1,1))/dy
     end if

     if(v(2,ny)==undeff.or.v(1,ny)==undeff.or.u(1,ny)==undeff.or.u(1,ny-1)==undeff)then
        val(1,ny)=undeff
     else
        val(1,ny)=(v(2,ny)-v(1,ny))/dx-(u(1,ny)-u(1,ny-1))/dy
     end if

     if(v(nx,1)==undeff.or.v(nx-1,1)==undeff.or.u(nx,2)==undeff.or.u(nx,1)==undeff)then
        val(nx,1)=undeff
     else
        val(nx,1)=(v(nx,1)-v(nx-1,1))/dx-(u(nx,2)-u(nx,1))/dy
     end if

     if(v(nx,ny)==undeff.or.v(nx-1,ny)==undeff.or.u(nx,ny)==undeff.or. u(nx,ny-1)==undeff)then
        val(nx,ny)=undeff
     else
        val(nx,ny)=(v(nx,ny)-v(nx-1,ny))/dx-(u(nx,ny)-u(nx,ny-1))/dy
     end if

  else

     do j=2,ny-1
        do i=2,nx-1
           val(i,j)=0.5*((v(i+1,j)-v(i-1,j))/dx-(u(i,j+1)-u(i,j-1))/dy)
        end do
     end do

!-- 領域の端 ---
!-- y 方向の両端 ---
     do i=2,nx-1
        val(i,1)=0.5*(v(i+1,1)-v(i-1,1))/dx-(u(i,2)-u(i,1))/dy
        val(i,ny)=0.5*(v(i+1,ny)-v(i-1,ny))/dx-(u(i,ny)-u(i,ny-1))/dy
     end do
!-- x 方向の両端 ---
     do j=2,ny-1
        val(1,j)=(v(2,j)-v(1,j))/dx-0.5*((u(1,j+1)-u(1,j-1))/dy)
        val(nx,j)=(v(nx,j)-v(nx-1,j))/dx-0.5*((u(nx,j+1)-u(nx,j-1))/dy)
     end do
!-- 4 隅 ---
     val(1,1)=(v(2,1)-v(1,1))/dx-(u(1,2)-u(1,1))/dy
     val(1,ny)=(v(2,ny)-v(1,ny))/dx-(u(1,ny)-u(1,ny-1))/dy
     val(nx,1)=(v(nx,1)-v(nx-1,1))/dx-(u(nx,2)-u(nx,1))/dy
     val(nx,ny)=(v(nx,ny)-v(nx-1,ny))/dx-(u(nx,ny)-u(nx,ny-1))/dy
  end if

  if(order==.false.)then
     do j=1,ny
        do i=1,nx
           val(i,j)=-val(i,j)
        end do
     end do
  end if

end subroutine curl
Subroutine :
x(:) :real, intent(in)
: x 方向の空間座標 [m]
y(:) :real, intent(in)
: y 方向の空間座標 [m]
u(size(x),size(y)) :real, intent(in)
: x に対応する方向の 2 次元ベクトル成分
v(size(x),size(y)) :real, intent(in)
: y に対応する方向の 2 次元ベクトル成分
val(size(x),size(y)) :real, intent(inout)
: 2次元発散値
undeff :real, intent(in), optional
hx(size(x),size(y)) :real, intent(in), optional
: x 方向のスケール因子
hy(size(x),size(y)) :real, intent(in), optional
: y 方向のスケール因子

2次元発散計算ルーチン

$frac{partial u}{partial x} +frac{partial v}{partial y} $ を 2 次の中央差分近似で書き換えると, 点 $(i,j)$ での発散は $frac{u_{i+1,j}-u_{i-1,j}}{2dx} + frac{v_{i,j+1}-v_{i,j-1}}{2dy} $ とできる. これを用いて2次元発散を計算. データ点が足りない隅の領域では, 1 次の差分近似で計算するので, 少し精度が 落ちる. 基本的に等間隔格子のデータのみしか計算できないので, 不等間隔のデータは 等間隔に内挿なりを行ってからこのルーチンを呼び出す必要あり. 将来, 任意の格子座標系に対応するため, 後ろにスケーリングファクターオプションを 追加したが, 現在は使用していない.

[Source]

subroutine div( x, y, u, v, val, undeff, hx, hy )
! 2次元発散計算ルーチン
! 
! $\frac{\partial u}{\partial x} +\frac{\partial v}{\partial y} $ を
! 2 次の中央差分近似で書き換えると, 点 $(i,j)$ での発散は
! $\frac{u_{i+1,j}-u_{i-1,j}}{2dx} + \frac{v_{i,j+1}-v_{i,j-1}}{2dy} $
! とできる. これを用いて2次元発散を計算.
! データ点が足りない隅の領域では, 1 次の差分近似で計算するので, 少し精度が
! 落ちる.
! 基本的に等間隔格子のデータのみしか計算できないので, 不等間隔のデータは
! 等間隔に内挿なりを行ってからこのルーチンを呼び出す必要あり.
! 将来, 任意の格子座標系に対応するため, 後ろにスケーリングファクターオプションを
! 追加したが, 現在は使用していない.
  implicit none
  real, intent(in) :: x(:)  ! x 方向の空間座標 [m]
  real, intent(in) :: y(:)  ! y 方向の空間座標 [m]
  real, intent(in) :: u(size(x),size(y))  ! x に対応する方向の 2 次元ベクトル成分
  real, intent(in) :: v(size(x),size(y))  ! y に対応する方向の 2 次元ベクトル成分
  real, intent(inout) :: val(size(x),size(y))  ! 2次元発散値
  real, intent(in), optional :: undeff
  real, intent(in), optional :: hx(size(x),size(y))  ! x 方向のスケール因子
  real, intent(in), optional :: hy(size(x),size(y))  ! y 方向のスケール因子
  integer :: i   ! イタレーション用添字
  integer :: j   ! イタレーション用添字
  integer :: nx  ! 空間配列要素数 1 次元目
  integer :: ny  ! 空間配列要素数 2 次元目
  real :: dx  ! 1 次元目を微分する格子間隔 [m]
  real :: dy  ! 2 次元目を微分する格子間隔 [m]

  dx=x(2)-x(1)
  dy=y(2)-y(1)
  nx=size(x)
  ny=size(y)

  if(present(undeff))then
     do j=2,ny-1
        do i=2,nx-1
           if(u(i+1,j)==undeff.or.u(i-1,j)==undeff.or.v(i,j+1)==undeff.or. v(i,j-1)==undeff)then
              val(i,j)=undeff
           else
              val(i,j)=0.5*((u(i+1,j)-u(i-1,j))/dx+(v(i,j+1)-v(i,j-1))/dy)
           end if
        end do
     end do

!-- 領域の端 ---
!-- y 方向の両端 ---
     do i=2,nx-1
        if(u(i+1,1)==undeff.or.u(i-1,1)==undeff.or.v(i,2)==undeff.or.v(i,1)==undeff)then
           val(i,1)=undeff
        else
           val(i,1)=0.5*(u(i+1,1)-u(i-1,1))/dx+(v(i,2)-v(i,1))/dy
        end if

        if(u(i+1,ny)==undeff.or.u(i-1,ny)==undeff.or.v(i,ny)==undeff.or. v(i,ny-1)==undeff)then
           val(i,ny)=undeff
        else
           val(i,ny)=0.5*(u(i+1,ny)-u(i-1,ny))/dx+(v(i,ny)-v(i,ny-1))/dy
        end if
     end do
!-- x 方向の両端 ---
     do j=2,ny-1
        if(u(2,j)==undeff.or.u(1,j)==undeff.or.v(1,j+1)==undeff.or.v(1,j-1)==undeff)then
           val(1,j)=undeff
        else
           val(1,j)=(u(2,j)-u(1,j))/dx+0.5*((v(1,j+1)-v(1,j-1))/dy)
        end if

        if(u(nx,j)==undeff.or.u(nx-1,j)==undeff.or.v(nx,j+1)==undeff.or. v(nx,j-1)==undeff)then
           val(nx,j)=undeff
        else
           val(nx,j)=(u(nx,j)-u(nx-1,j))/dx+0.5*((v(nx,j+1)-v(nx,j-1))/dy)
        end if
     end do
!-- 4 隅 ---
     if(u(2,1)==undeff.or.u(1,1)==undeff.or.v(1,2)==undeff.or.v(1,1)==undeff)then
        val(1,1)=undeff
     else
        val(1,1)=(u(2,1)-u(1,1))/dx+(v(1,2)-v(1,1))/dy
     end if

     if(u(2,ny)==undeff.or.u(1,ny)==undeff.or.v(1,ny)==undeff.or.v(1,ny-1)==undeff)then
        val(1,ny)=undeff
     else
        val(1,ny)=(u(2,ny)-u(1,ny))/dx+(v(1,ny)-v(1,ny-1))/dy
     end if

     if(u(nx,1)==undeff.or.u(nx-1,1)==undeff.or.v(nx,2)==undeff.or.v(nx,1)==undeff)then
        val(nx,1)=undeff
     else
        val(nx,1)=(u(nx,1)-u(nx-1,1))/dx+(v(nx,2)-v(nx,1))/dy
     end if

     if(u(nx,ny)==undeff.or.u(nx-1,ny)==undeff.or.v(nx,ny)==undeff.or. v(nx,ny-1)==undeff)then
        val(nx,ny)=undeff
     else
        val(nx,ny)=(u(nx,ny)-u(nx-1,ny))/dx+(v(nx,ny)-v(nx,ny-1))/dy
     end if

  else

     do j=2,ny-1
        do i=2,nx-1
           val(i,j)=0.5*((u(i+1,j)-u(i-1,j))/dx+(v(i,j+1)-v(i,j-1))/dy)
        end do
     end do

!-- 領域の端 ---
!-- y 方向の両端 ---
     do i=2,nx-1
        val(i,1)=0.5*(u(i+1,1)-u(i-1,1))/dx+(v(i,2)-v(i,1))/dy
        val(i,ny)=0.5*(u(i+1,ny)-u(i-1,ny))/dx+(v(i,ny)-v(i,ny-1))/dy
     end do
!-- x 方向の両端 ---
     do j=2,ny-1
        val(1,j)=(u(2,j)-u(1,j))/dx+0.5*((v(1,j+1)-v(1,j-1))/dy)
        val(nx,j)=(u(nx,j)-u(nx-1,j))/dx+0.5*((v(nx,j+1)-v(nx,j-1))/dy)
     end do
!-- 4 隅 ---
     val(1,1)=(u(2,1)-u(1,1))/dx+(v(1,2)-v(1,1))/dy
     val(1,ny)=(u(2,ny)-u(1,ny))/dx+(v(1,ny)-v(1,ny-1))/dy
     val(nx,1)=(u(nx,1)-u(nx-1,1))/dx+(v(nx,2)-v(nx,1))/dy
     val(nx,ny)=(u(nx,ny)-u(nx-1,ny))/dx+(v(nx,ny)-v(nx,ny-1))/dy

  end if

end subroutine div
Subroutine :
x(:) :real, intent(in)
: x 方向の空間座標 [m]
y(:) :real, intent(in)
: y 方向の空間座標 [m]
z(:) :real, intent(in)
: z 方向の空間座標 [m]
u(size(x),size(y),size(z)) :real, intent(in)
: x に対応する方向の 2 次元ベクトル成分
v(size(x),size(y),size(z)) :real, intent(in)
: y に対応する方向の 2 次元ベクトル成分
w(size(x),size(y),size(z)) :real, intent(in)
: y に対応する方向の 2 次元ベクトル成分
val(size(x),size(y),size(z)) :real, intent(inout)
: 3 次元発散値
undeff :real, intent(in), optional
hx(size(x),size(y),size(z)) :real, intent(in), optional
: x 方向のスケール因子
hy(size(x),size(y),size(z)) :real, intent(in), optional
: y 方向のスケール因子
hz(size(x),size(y),size(z)) :real, intent(in), optional
: z 方向のスケール因子

3次元発散計算ルーチン

$frac{partial u}{partial x} +frac{partial v}{partial y} +frac{partial w}{partial z} $ を 2 次の中央差分近似で書き換えると, 点 $(i,j,k)$ での発散は $frac{u_{i+1,j,k}-u_{i-1,j,k}}{2dx} + frac{v_{i,j+1,k}-v_{i,j-1,k}}{2dy} + frac{w_{i,j,k+1}-w_{i,j,k-1}}{2dz} $ とできる. これを用いて2次元発散を計算. データ点が足りない隅の領域では, 1 次の差分近似で計算するので, 少し精度が 落ちる. 基本的に等間隔格子のデータのみしか計算できないので, 不等間隔のデータは 等間隔に内挿なりを行ってからこのルーチンを呼び出す必要あり. 将来, 任意の格子座標系に対応するため, 後ろにスケーリングファクターオプションを 追加したが, 現在は使用していない.

[Source]

subroutine div_3d( x, y, z, u, v, w, val, undeff, hx, hy, hz )
! 3次元発散計算ルーチン
!
! $\frac{\partial u}{\partial x} +\frac{\partial v}{\partial y} +\frac{\partial w}{\partial z} $ を
! 2 次の中央差分近似で書き換えると, 点 $(i,j,k)$ での発散は
! $\frac{u_{i+1,j,k}-u_{i-1,j,k}}{2dx} + \frac{v_{i,j+1,k}-v_{i,j-1,k}}{2dy} + \frac{w_{i,j,k+1}-w_{i,j,k-1}}{2dz} $
! とできる. これを用いて2次元発散を計算.
! データ点が足りない隅の領域では, 1 次の差分近似で計算するので, 少し精度が
! 落ちる.
! 基本的に等間隔格子のデータのみしか計算できないので, 不等間隔のデータは
! 等間隔に内挿なりを行ってからこのルーチンを呼び出す必要あり.
! 将来, 任意の格子座標系に対応するため, 後ろにスケーリングファクターオプションを
! 追加したが, 現在は使用していない.
  implicit none
  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))  ! x に対応する方向の 2 次元ベクトル成分
  real, intent(in) :: v(size(x),size(y),size(z))  ! y に対応する方向の 2 次元ベクトル成分
  real, intent(in) :: w(size(x),size(y),size(z))  ! y に対応する方向の 2 次元ベクトル成分
  real, intent(inout) :: val(size(x),size(y),size(z))  ! 3 次元発散値
  real, intent(in), optional :: undeff
  real, intent(in), optional :: hx(size(x),size(y),size(z))  ! x 方向のスケール因子
  real, intent(in), optional :: hy(size(x),size(y),size(z))  ! y 方向のスケール因子
  real, intent(in), optional :: hz(size(x),size(y),size(z))  ! z 方向のスケール因子
  integer :: i   ! イタレーション用添字
  integer :: j   ! イタレーション用添字
  integer :: k   ! イタレーション用添字
  integer :: nx  ! 空間配列要素数 1 次元目
  integer :: ny  ! 空間配列要素数 2 次元目
  integer :: nz  ! 空間配列要素数 3 次元目
  real :: dx  ! 1 次元目を微分する格子間隔 [m]
  real :: dy  ! 2 次元目を微分する格子間隔 [m]
  real :: dz  ! 3 次元目を微分する格子間隔 [m]

  dx=x(2)-x(1)
  dy=y(2)-y(1)
  dz=z(2)-z(1)
  nx=size(x)
  ny=size(y)
  nz=size(z)

  if(present(undeff))then
!$omp parallel default(shared)
!$omp do private(i,j,k)
     do k=2,nz-1
        do j=2,ny-1
           do i=2,nx-1
              if(u(i+1,j,k)==undeff.or.u(i-1,j,k)==undeff.or. v(i,j+1,k)==undeff.or.v(i,j-1,k)==undeff.or. w(i,j,k+1)==undeff.or.w(i,j,k-1)==undeff)then
                 val(i,j,k)=undeff
              else
                 val(i,j,k)=0.5*((u(i+1,j,k)-u(i-1,j,k))/dx +(v(i,j+1,k)-v(i,j-1,k))/dy +(w(i,j,k+1)-w(i,j,k-1))/dz)
              end if
           end do
        end do
     end do
!$omp end do
!$omp end parallel
  else
!$omp parallel default(shared)
!$omp do private(i,j,k)
     do k=2,nz-1
        do j=2,ny-1
           do i=2,nx-1
              val(i,j,k)=0.5*((u(i+1,j,k)-u(i-1,j,k))/dx +(v(i,j+1,k)-v(i,j-1,k))/dy +(w(i,j,k+1)-w(i,j,k-1))/dz)
           end do
        end do
     end do
!$omp end do
!$omp end parallel
  end if

!-- 領域の端を除く面 ---
!-- この面の計算は div ルーチンで先に計算させ, その結果に面に垂直な成分の
!   発散を足し合わせることで計算させる.
!-- z 軸に垂直な面 ---
  if(present(undeff))then
     call div( x, y, u(:,:,1), v(:,:,1), val(:,:,1), undeff )
     call div( x, y, u(:,:,nz), v(:,:,nz), val(:,:,nz), undeff )
     do j=2,ny-1
        do i=2,nx-1
           if(val(i,j,1)==undeff.or.w(i,j,2)==undeff.or.w(i,j,1)==undeff)then
              val(i,j,1)=undeff
           else
              val(i,j,1)=val(i,j,1)+(w(i,j,2)-w(i,j,1))/dz
           end if

           if(val(i,j,nz)==undeff.or.w(i,j,nz)==undeff.or.w(i,j,nz-1)==undeff)then
              val(i,j,nz)=undeff
           else
              val(i,j,nz)=val(i,j,nz)+(w(i,j,nz)-w(i,j,nz-1))/dz
           end if
        end do
     end do
  else
     call div( x, y, u(:,:,1), v(:,:,1), val(:,:,1) )
     call div( x, y, u(:,:,nz), v(:,:,nz), val(:,:,nz) )
     do j=2,ny-1
        do i=2,nx-1
           val(i,j,1)=val(i,j,1)+(w(i,j,2)-w(i,j,1))/dz
           val(i,j,nz)=val(i,j,nz)+(w(i,j,nz)-w(i,j,nz-1))/dz
        end do
     end do
  end if

!-- y 軸に垂直な面 ---
  if(present(undeff))then
     call div( x, z, u(:,1,:), w(:,1,:), val(:,1,:), undeff )
     call div( x, z, u(:,ny,:), w(:,ny,:), val(:,ny,:), undeff )
     do k=2,nz-1
        do i=2,nx-1
           if(val(i,1,k)==undeff.or.v(i,2,k)==undeff.or.v(i,1,k)==undeff)then
              val(i,1,k)=undeff
           else
              val(i,1,k)=val(i,1,k)+(v(i,2,k)-v(i,1,k))/dy
           end if

           if(val(i,ny,k)==undeff.or.v(i,ny,k)==undeff.or.v(i,ny-1,k)==undeff)then
              val(i,ny,k)=undeff
           else
              val(i,ny,k)=val(i,ny,k)+(v(i,ny,k)-v(i,ny-1,k))/dy
           end if
        end do
     end do
  else
     call div( x, z, u(:,1,:), w(:,1,:), val(:,1,:) )
     call div( x, z, u(:,ny,:), w(:,ny,:), val(:,ny,:) )
     do k=2,nz-1
        do i=2,nx-1
           val(i,1,k)=val(i,1,k)+(v(i,2,k)-v(i,1,k))/dy
           val(i,ny,k)=val(i,ny,k)+(v(i,ny,k)-v(i,ny-1,k))/dy
        end do
     end do
  end if

!-- x 軸に垂直な面 ---
  if(present(undeff))then
     call div( y, z, v(1,:,:), w(1,:,:), val(1,:,:), undeff )
     call div( y, z, v(nx,:,:), w(nx,:,:), val(nx,:,:), undeff )
     do k=2,nz-1
        do j=2,ny-1
           if(val(1,j,k)==undeff.or.u(2,j,k)==undeff.or.u(1,j,k)==undeff)then
              val(1,j,k)=undeff
           else
              val(1,j,k)=val(1,j,k)+(u(2,j,k)-u(1,j,k))/dx
           end if

           if(val(nx,j,k)==undeff.or.u(nx,j,k)==undeff.or.u(nx-1,j,k)==undeff)then
              val(nx,j,k)=undeff
           else
              val(nx,j,k)=val(nx,j,k)+(u(nx,j,k)-u(nx-1,j,k))/dx
           end if
        end do
     end do
  else
     call div( y, z, v(1,:,:), w(1,:,:), val(1,:,:) )
     call div( y, z, v(nx,:,:), w(nx,:,:), val(nx,:,:) )
     do k=2,nz-1
        do j=2,ny-1
           val(1,j,k)=val(1,j,k)+(u(2,j,k)-u(1,j,k))/dx
           val(nx,j,k)=val(nx,j,k)+(u(nx,j,k)-u(nx-1,j,k))/dx
        end do
     end do
  end if

!-- 各軸に平行な境界線上 12 本
!-- x 軸に平行な 4 本
  if(present(undeff))then
     do i=2,nx-1
        if(u(i+1,1,1)==undeff.or.u(i-1,1,1)==undeff.or. v(i,2,1)==undeff.or.v(i,1,1)==undeff.or. w(i,1,2)==undeff.or.w(i,1,1)==undeff)then
           val(i,1,1)=undeff
        else
           val(i,1,1)=0.5*(u(i+1,1,1)-u(i-1,1,1))/dx +(v(i,2,1)-v(i,1,1))/dy +(w(i,1,2)-w(i,1,1))/dz
        end if

        if(u(i+1,ny,1)==undeff.or.u(i-1,ny,1)==undeff.or. v(i,ny,1)==undeff.or.v(i,ny-1,1)==undeff.or. w(i,ny,2)==undeff.or.w(i,ny,1)==undeff)then
           val(i,ny,1)=undeff
        else
           val(i,ny,1)=0.5*(u(i+1,ny,1)-u(i-1,ny,1))/dx +(v(i,ny,1)-v(i,ny-1,1))/dy +(w(i,ny,2)-w(i,ny,1))/dz
        end if

        if(u(i+1,1,nz)==undeff.or.u(i-1,1,nz)==undeff.or. v(i,2,nz)==undeff.or.v(i,1,nz)==undeff.or. w(i,1,nz)==undeff.or.w(i,1,nz-1)==undeff)then
           val(i,1,nz)=undeff
        else
           val(i,1,nz)=0.5*(u(i+1,1,nz)-u(i-1,1,nz))/dx +(v(i,2,nz)-v(i,1,nz))/dy +(w(i,1,nz)-w(i,1,nz-1))/dz
        end if

        if(u(i+1,ny,nz)==undeff.or.u(i-1,ny,nz)==undeff.or. v(i,ny,nz)==undeff.or.v(i,ny-1,nz)==undeff.or. w(i,ny,nz)==undeff.or.w(i,ny,nz-1)==undeff)then
           val(i,ny,nz)=undeff
        else
           val(i,ny,nz)=0.5*(u(i+1,ny,nz)-u(i-1,ny,nz))/dx +(v(i,ny,nz)-v(i,ny-1,nz))/dy +(w(i,ny,nz)-w(i,ny,nz-1))/dz
        end if
     end do
  else
     val(i,1,1)=0.5*(u(i+1,1,1)-u(i-1,1,1))/dx +(v(i,2,1)-v(i,1,1))/dy +(w(i,1,2)-w(i,1,1))/dz
     val(i,ny,1)=0.5*(u(i+1,ny,1)-u(i-1,ny,1))/dx +(v(i,ny,1)-v(i,ny-1,1))/dy +(w(i,ny,2)-w(i,ny,1))/dz
     val(i,1,nz)=0.5*(u(i+1,1,nz)-u(i-1,1,nz))/dx +(v(i,2,nz)-v(i,1,nz))/dy +(w(i,1,nz)-w(i,1,nz-1))/dz
     val(i,ny,nz)=0.5*(u(i+1,ny,nz)-u(i-1,ny,nz))/dx +(v(i,ny,nz)-v(i,ny-1,nz))/dy +(w(i,ny,nz)-w(i,ny,nz-1))/dz
  end if

!-- y 軸に平行な 4 本
  if(present(undeff))then
     do j=2,ny-1
        if(u(2,j,1)==undeff.or.u(1,j,1)==undeff.or. v(1,j+1,1)==undeff.or.v(1,j-1,1)==undeff.or. w(1,j,2)==undeff.or.w(1,j,1)==undeff)then
           val(1,j,1)=undeff
        else
           val(1,j,1)=(u(2,j,1)-u(1,j,1))/dx +0.5*(v(1,j+1,1)-v(1,j-1,1))/dy +(w(1,j,2)-w(1,j,1))/dz
        end if

        if(u(nx,j,1)==undeff.or.u(nx-1,j,1)==undeff.or. v(nx,j+1,1)==undeff.or.v(nx,j-1,1)==undeff.or. w(nx,j,2)==undeff.or.w(nx,j,1)==undeff)then
           val(nx,j,1)=undeff
        else
           val(nx,j,1)=(u(nx,j,1)-u(nx-1,j,1))/dx +0.5*(v(nx,j+1,1)-v(nx,j-1,1))/dy +(w(nx,j,2)-w(nx,j,1))/dz
        end if

        if(u(2,j,nz)==undeff.or.u(1,j,nz)==undeff.or. v(1,j+1,nz)==undeff.or.v(1,j-1,nz)==undeff.or. w(1,j,nz)==undeff.or.w(1,j,nz-1)==undeff)then
           val(1,j,nz)=undeff
        else
           val(1,j,nz)=(u(2,j,nz)-u(1,j,nz))/dx +0.5*(v(1,j+1,nz)-v(1,j-1,nz))/dy +(w(1,j,nz)-w(1,j,nz-1))/dz
        end if

        if(u(nx,j,nz)==undeff.or.u(nx-1,j,nz)==undeff.or. v(nx,j+1,nz)==undeff.or.v(nx,j-1,nz)==undeff.or. w(nx,j,nz)==undeff.or.w(nx,j,nz-1)==undeff)then
           val(nx,j,nz)=undeff
        else
           val(nx,j,nz)=(u(nx,j,nz)-u(nx-1,j,nz))/dx +0.5*(v(nx,j+1,nz)-v(nx,j-1,nz))/dy +(w(nx,j,nz)-w(nx,j,nz-1))/dz
        end if
     end do
  else
     do j=2,ny-1
        val(1,j,1)=(u(2,j,1)-u(1,j,1))/dx +0.5*(v(1,j+1,1)-v(1,j-1,1))/dy +(w(1,j,2)-w(1,j,1))/dz
        val(nx,j,1)=(u(nx,j,1)-u(nx-1,j,1))/dx +0.5*(v(nx,j+1,1)-v(nx,j-1,1))/dy +(w(nx,j,2)-w(nx,j,1))/dz
        val(1,j,nz)=(u(2,j,nz)-u(1,j,nz))/dx +0.5*(v(1,j+1,nz)-v(1,j-1,nz))/dy +(w(1,j,nz)-w(1,j,nz-1))/dz
        val(nx,j,nz)=(u(nx,j,nz)-u(nx-1,j,nz))/dx +0.5*(v(nx,j+1,nz)-v(nx,j-1,nz))/dy +(w(nx,j,nz)-w(nx,j,nz-1))/dz
     end do
  end if

!-- z 軸に平行な 4 本
  if(present(undeff))then
     do k=2,nz-1
        if(u(2,1,k)==undeff.or.u(1,1,k)==undeff.or. v(1,2,k)==undeff.or.v(1,1,k)==undeff.or. w(1,1,k+1)==undeff.or.w(1,1,k-1)==undeff)then
           val(1,1,k)=undeff
        else
           val(1,1,k)=(u(2,1,k)-u(1,1,k))/dx +(v(1,2,k)-v(1,1,k))/dy +0.5*(w(1,1,k+1)-w(1,1,k-1))/dz
        end if

        if(u(nx,1,k)==undeff.or.u(nx-1,1,k)==undeff.or. v(nx,2,k)==undeff.or.v(nx,1,k)==undeff.or. w(nx,1,k+1)==undeff.or.w(nx,1,k-1)==undeff)then
           val(nx,1,k)=undeff
        else
           val(nx,1,k)=(u(nx,1,k)-u(nx-1,1,k))/dx +(v(nx,2,k)-v(nx,1,k))/dy +0.5*(w(nx,1,k+1)-w(nx,1,k-1))/dz
        end if

        if(u(2,ny,k)==undeff.or.u(1,ny,k)==undeff.or. v(1,ny,k)==undeff.or.v(1,ny-1,k)==undeff.or. w(1,ny,k+1)==undeff.or.w(1,ny,k-1)==undeff)then
           val(1,ny,k)=undeff
        else
           val(1,ny,k)=(u(2,ny,k)-u(1,ny,k))/dx +(v(1,ny,k)-v(1,ny-1,k))/dy +0.5*(w(1,ny,k+1)-w(1,ny,k-1))/dz
        end if

        if(u(nx,ny,k)==undeff.or.u(nx-1,ny,k)==undeff.or. v(nx,ny,k)==undeff.or.v(nx,ny-1,k)==undeff.or. w(nx,ny,k+1)==undeff.or.w(nx,ny,k-1)==undeff)then
           val(nx,ny,k)=undeff
        else
           val(nx,ny,k)=(u(nx,ny,k)-u(nx-1,ny,k))/dx +(v(nx,ny,k)-v(nx,ny-1,k))/dy +0.5*(w(nx,ny,k+1)-w(nx,ny,k-1))/dz
        end if
     end do
  else
     do k=2,nz-1
        val(1,1,k)=(u(2,1,k)-u(1,1,k))/dx +(v(1,2,k)-v(1,1,k))/dy +0.5*(w(1,1,k+1)-w(1,1,k-1))/dz
        val(nx,1,k)=(u(nx,1,k)-u(nx-1,1,k))/dx +(v(nx,2,k)-v(nx,1,k))/dy +0.5*(w(nx,1,k+1)-w(nx,1,k-1))/dz
        val(1,ny,k)=(u(2,ny,k)-u(1,ny,k))/dx +(v(1,ny,k)-v(1,ny-1,k))/dy +0.5*(w(1,ny,k+1)-w(1,ny,k-1))/dz
        val(nx,ny,k)=(u(nx,ny,k)-u(nx-1,ny,k))/dx +(v(nx,ny,k)-v(nx,ny-1,k))/dy +0.5*(w(nx,ny,k+1)-w(nx,ny,k-1))/dz
     end do
  end if

!-- 8 隅 ---
  if(present(undeff))then
     if(u(2,1,1)==undeff.or.-u(1,1,1)==undeff.or. v(1,2,1)==undeff.or.v(1,1,1)==undeff.or. w(1,1,2)==undeff.or.w(1,1,1)==undeff)then
        val(1,1,1)=undeff
     else
        val(1,1,1)=(u(2,1,1)-u(1,1,1))/dx +(v(1,2,1)-v(1,1,1))/dy +(w(1,1,2)-w(1,1,1))/dz
     end if

     if(u(nx,1,1)==undeff.or.u(nx-1,1,1)==undeff.or. v(nx,2,1)==undeff.or.v(nx,1,1)==undeff.or. w(nx,1,2)==undeff.or.w(nx,1,1)==undeff)then
        val(nx,1,1)=undeff
     else
        val(nx,1,1)=(u(nx,1,1)-u(nx-1,1,1))/dx +(v(nx,2,1)-v(nx,1,1))/dy +(w(nx,1,2)-w(nx,1,1))/dz
     end if

     if(u(2,ny,1)==undeff.or.u(1,ny,1)==undeff.or. v(1,ny,1)==undeff.or.v(1,ny-1,1)==undeff.or. w(1,ny,2)==undeff.or.w(1,ny,1)==undeff)then
        val(1,ny,1)=undeff
     else
        val(1,ny,1)=(u(2,ny,1)-u(1,ny,1))/dx +(v(1,ny,1)-v(1,ny-1,1))/dy +(w(1,ny,2)-w(1,ny,1))/dz
     end if

     if(u(2,1,nz)==undeff.or.u(1,1,nz)==undeff.or. v(1,2,nz)==undeff.or.v(1,1,nz)==undeff.or. w(1,1,nz)==undeff.or.w(1,1,nz-1)==undeff)then
        val(1,1,nz)=undeff
     else
        val(1,1,nz)=(u(2,1,nz)-u(1,1,nz))/dx +(v(1,2,nz)-v(1,1,nz))/dy +(w(1,1,nz)-w(1,1,nz-1))/dz
     end if

     if(u(nx,ny,1)==undeff.or.u(nx-1,ny,1)==undeff.or. v(nx,ny,1)==undeff.or.v(nx,ny-1,1)==undeff.or. w(nx,ny,2)==undeff.or.w(nx,ny,1)==undeff)then
        val(nx,ny,1)=undeff
     else
        val(nx,ny,1)=(u(nx,ny,1)-u(nx-1,ny,1))/dx +(v(nx,ny,1)-v(nx,ny-1,1))/dy +(w(nx,ny,2)-w(nx,ny,1))/dz
     end if

     if(u(2,ny,nz)==undeff.or.u(1,ny,nz)==undeff.or. v(1,ny,nz)==undeff.or.v(1,ny-1,nz)==undeff.or. w(1,ny,nz)==undeff.or.w(1,ny,nz-1)==undeff)then
        val(1,ny,nz)=undeff
     else
        val(1,ny,nz)=(u(2,ny,nz)-u(1,ny,nz))/dx +(v(1,ny,nz)-v(1,ny-1,nz))/dy +(w(1,ny,nz)-w(1,ny,nz-1))/dz
     end if

     if(u(nx,1,nz)==undeff.or.u(nx-1,1,nz)==undeff.or. v(nx,2,nz)==undeff.or.v(nx,1,nz)==undeff.or. w(nx,1,nz)==undeff.or.w(nx,1,nz-1)==undeff)then
        val(nx,1,nz)=undeff
     else
        val(nx,1,nz)=(u(nx,1,nz)-u(nx-1,1,nz))/dx +(v(nx,2,nz)-v(nx,1,nz))/dy +(w(nx,1,nz)-w(nx,1,nz-1))/dz
     end if

     if(u(nx,ny,nz)==undeff.or.u(nx-1,ny,nz)==undeff.or. v(nx,ny,nz)==undeff.or.v(nx,ny-1,nz)==undeff.or. w(nx,ny,nz)==undeff.or.w(nx,ny,nz-1)==undeff)then
        val(nx,ny,nz)=undeff
     else
        val(nx,ny,nz)=(u(nx,ny,nz)-u(nx-1,ny,nz))/dx +(v(nx,ny,nz)-v(nx,ny-1,nz))/dy +(w(nx,ny,nz)-w(nx,ny,nz-1))/dz
     end if
  else
     val(1,1,1)=(u(2,1,1)-u(1,1,1))/dx +(v(1,2,1)-v(1,1,1))/dy +(w(1,1,2)-w(1,1,1))/dz
     val(nx,1,1)=(u(nx,1,1)-u(nx-1,1,1))/dx +(v(nx,2,1)-v(nx,1,1))/dy +(w(nx,1,2)-w(nx,1,1))/dz
     val(1,ny,1)=(u(2,ny,1)-u(1,ny,1))/dx +(v(1,ny,1)-v(1,ny-1,1))/dy +(w(1,ny,2)-w(1,ny,1))/dz
     val(1,1,nz)=(u(2,1,nz)-u(1,1,nz))/dx +(v(1,2,nz)-v(1,1,nz))/dy +(w(1,1,nz)-w(1,1,nz-1))/dz
     val(nx,ny,1)=(u(nx,ny,1)-u(nx-1,ny,1))/dx +(v(nx,ny,1)-v(nx,ny-1,1))/dy +(w(nx,ny,2)-w(nx,ny,1))/dz
     val(1,ny,nz)=(u(2,ny,nz)-u(1,ny,nz))/dx +(v(1,ny,nz)-v(1,ny-1,nz))/dy +(w(1,ny,nz)-w(1,ny,nz-1))/dz
     val(nx,1,nz)=(u(nx,1,nz)-u(nx-1,1,nz))/dx +(v(nx,2,nz)-v(nx,1,nz))/dy +(w(nx,1,nz)-w(nx,1,nz-1))/dz
     val(nx,ny,nz)=(u(nx,ny,nz)-u(nx-1,ny,nz))/dx +(v(nx,ny,nz)-v(nx,ny-1,nz))/dy +(w(nx,ny,nz)-w(nx,ny,nz-1))/dz
  end if

end subroutine
Subroutine :
x(:,:,:) :real, intent(in)
: x 方向のベクトル成分
y(size(x,1),size(x,2),size(x,3)) :real, intent(in)
: y 方向のベクトル成分
z(size(x,1),size(x,2),size(x,3)) :real, intent(in)
: z 方向のベクトル成分
u(size(x,1),size(x,2),size(x,3)) :real, intent(in)
: x 方向のベクトル成分
v(size(x,1),size(x,2),size(x,3)) :real, intent(in)
: y 方向のベクトル成分
w(size(x,1),size(x,2),size(x,3)) :real, intent(in)
: z 方向のベクトル成分
dot(size(x,1),size(x,2),size(x,3)) :real, intent(inout)
: 内積
undeff :real, intent(in), optional

2ベクトルの内積計算ルーチン 配列を工夫すると, 2 次元での内積を計算することも可能

[Source]

subroutine dot_prod(x,y,z,u,v,w,dot,undeff)
  ! 2ベクトルの内積計算ルーチン
  ! 配列を工夫すると, 2 次元での内積を計算することも可能
  implicit none
  real, intent(in) :: x(:,:,:)  ! x 方向のベクトル成分
  real, intent(in) :: y(size(x,1),size(x,2),size(x,3))  ! y 方向のベクトル成分
  real, intent(in) :: z(size(x,1),size(x,2),size(x,3))  ! z 方向のベクトル成分
  real, intent(in) :: u(size(x,1),size(x,2),size(x,3))  ! x 方向のベクトル成分
  real, intent(in) :: v(size(x,1),size(x,2),size(x,3))  ! y 方向のベクトル成分
  real, intent(in) :: w(size(x,1),size(x,2),size(x,3))  ! z 方向のベクトル成分
  real, intent(inout) :: dot(size(x,1),size(x,2),size(x,3))  ! 内積
  real, intent(in), optional :: undeff
  integer :: i, j, k, nx, ny, nz

  nx=size(x,1)
  ny=size(x,2)
  nz=size(x,3)

  if(present(undeff))then
!$omp parallel do shared(dot,x,y,z,u,v,w) private(i,j,k)
     do k=1,nz
        do j=1,ny
           do i=1,nx
              if(x(i,j,k)==undeff.or.u(i,j,k)==undeff.or.y(i,j,k)==undeff.or. v(i,j,k)==undeff.or.z(i,j,k)==undeff.or.w(i,j,k)==undeff)then
                 dot(i,j,k)=undeff
              else
                 dot(i,j,k)=x(i,j,k)*u(i,j,k)+y(i,j,k)*v(i,j,k)+z(i,j,k)*w(i,j,k)
              end if
           end do
        end do
     end do
!$omp end parallel do
  else
!$omp parallel do shared(dot,x,y,z,u,v,w) private(i,j,k)
     do k=1,nz
        do j=1,ny
           do i=1,nx
              dot(i,j,k)=x(i,j,k)*u(i,j,k)+y(i,j,k)*v(i,j,k)+z(i,j,k)*w(i,j,k)
           end do
        end do
     end do
!$omp end parallel do
  end if

end subroutine dot_prod
Subroutine :
x(:) :real, intent(in)
: 空間座標
u(size(x)) :real, intent(in)
: 上の空間配列に対応する 1 次元スカラー値
val(size(x)) :real, intent(inout)
: スカラー値の x 方向の勾配
undef :real, intent(in), optional
hx(size(x)) :real, intent(in), optional
: x のスケール因子

1 次元のスカラー変数の勾配を計算する $frac{partial p}{partial x} $ を 2 次の中央差分近似で書き換えると, 点 $(i)$ での勾配は $frac{p_{i+1}-p_{i-1}}{2dx} $ とできる. これを用いて 1 次元勾配を計算. データ点が足りない隅の領域では, 1 次の差分近似で計算するので, 少し精度が 落ちる. 基本的に等間隔格子のデータのみしか計算できないので, 不等間隔のデータは 等間隔に内挿なりを行ってからこのルーチンを呼び出す必要あり.

[Source]

subroutine grad_1d( x, u, val, undef, hx )
! 1 次元のスカラー変数の勾配を計算する
! $\frac{\partial p}{\partial x} $ を
! 2 次の中央差分近似で書き換えると, 点 $(i)$ での勾配は
! $\frac{p_{i+1}-p_{i-1}}{2dx} $
! とできる. これを用いて 1 次元勾配を計算.
! データ点が足りない隅の領域では, 1 次の差分近似で計算するので, 少し精度が
! 落ちる.
! 基本的に等間隔格子のデータのみしか計算できないので, 不等間隔のデータは
! 等間隔に内挿なりを行ってからこのルーチンを呼び出す必要あり.
  implicit none
  real, intent(in) :: x(:)  ! 空間座標
  real, intent(in) :: u(size(x))  ! 上の空間配列に対応する 1 次元スカラー値
  real, intent(inout) :: val(size(x))  ! スカラー値の x 方向の勾配
  real, intent(in), optional :: undef
  real, intent(in), optional :: hx(size(x))  ! x のスケール因子
  integer :: i  ! イタレーション用添字
  integer :: nx  ! 配列要素数
  real :: dx  ! 格子間隔 [m]

  nx=size(x)
  dx=x(2)-x(1)

  if(present(undef))then
     do i=2,nx-1
        if(u(i+1)==undef.or.u(i-1)==undef)then
           val(i)=undef
        else
           val(i)=0.5*(u(i+1)-u(i-1))/dx
        end if
     end do
!-- データ数のない両端の処理 ---
     if(u(1)==undef.or.u(2)==undef)then
        val(1)=undef
     else
        val(1)=(u(2)-u(1))/dx
     end if
     if(u(nx)==undef.or.u(nx-1)==undef)then
        val(nx)=undef
     else
        val(nx)=(u(nx)-u(nx-1))/dx
     end if
  else
     do i=2,nx-1
        val(i)=0.5*(u(i+1)-u(i-1))/dx
     end do
!-- データ数のない両端の処理 ---
     val(1)=(u(2)-u(1))/dx
     val(nx)=(u(nx)-u(nx-1))/dx
  end if

end subroutine grad_1d
Subroutine :
x(:) :real, intent(in)
: x 方向の座標変数 [m]
y(:) :real, intent(in)
: y 方向の座標変数 [m]
u(size(x),size(y)) :real, intent(in)
: 勾配をとる 2 次元スカラー成分
valx(size(x),size(y)) :real, intent(inout)
: 計算された y 方向の 2 次元勾配ベクトル
valy(size(x),size(y)) :real, intent(inout)
: 計算された y 方向の 2 次元勾配ベクトル
undeff :real, intent(in), optional
hx(size(x),size(y)) :real, intent(in), optional
: x 方向のスケール因子
hy(size(x),size(y)) :real, intent(in), optional
: y 方向のスケール因子

$nabla _hp =left(frac{partial p}{partial x} ,\; frac{partial p}{partial y} right) $ を 2 次の中央差分近似で書き換えると, 点 $(i,j)$ での勾配は $left(frac{p_{i+1,j}-p_{i-1,j}}{2dx} ,\; frac{p_{i,j+1}-p_{i,j-1}}{2dy} right) $ とできる. これを用いて2次元勾配を計算. データ点が足りない隅の領域では, 1 次の差分近似で計算するので, 少し精度が 落ちる. 基本的に等間隔格子のデータのみしか計算できないので, 不等間隔のデータは 等間隔に内挿なりを行ってからこのルーチンを呼び出す必要あり. 先に用いた 1 次元勾配計算ルーチンを 2 回呼び出すことにしている.


[Source]

subroutine grad_vec_2d( x, y, u, valx, valy, undeff, hx, hy )
  !-------------------------------------------------------------------
  ! $\nabla _hp =\left(\frac{\partial p}{\partial x} ,\; \frac{\partial p}{\partial y} \right) $ を
  ! 2 次の中央差分近似で書き換えると, 点 $(i,j)$ での勾配は
  ! $\left(\frac{p_{i+1,j}-p_{i-1,j}}{2dx} ,\; \frac{p_{i,j+1}-p_{i,j-1}}{2dy} \right) $
  ! とできる. これを用いて2次元勾配を計算.
  ! データ点が足りない隅の領域では, 1 次の差分近似で計算するので, 少し精度が
  ! 落ちる.
  ! 基本的に等間隔格子のデータのみしか計算できないので, 不等間隔のデータは
  ! 等間隔に内挿なりを行ってからこのルーチンを呼び出す必要あり.
  ! 先に用いた 1 次元勾配計算ルーチンを 2 回呼び出すことにしている.
  !-------------------------------------------------------------------
  implicit none
  real, intent(in) :: x(:)  ! x 方向の座標変数 [m]
  real, intent(in) :: y(:)  ! y 方向の座標変数 [m]
  real, intent(in) :: u(size(x),size(y))  ! 勾配をとる 2 次元スカラー成分
  real, intent(inout) :: valx(size(x),size(y))  ! 計算された y 方向の 2 次元勾配ベクトル
  real, intent(inout) :: valy(size(x),size(y))  ! 計算された y 方向の 2 次元勾配ベクトル
  real, intent(in), optional :: undeff
  real, intent(in), optional :: hx(size(x),size(y))  ! x 方向のスケール因子
  real, intent(in), optional :: hy(size(x),size(y))  ! y 方向のスケール因子
  integer :: i   ! イタレーション用添字
  integer :: j   ! イタレーション用添字
  integer :: nx  ! x の配列要素数(size 関数で自動的に計算)
  integer :: ny  ! y の配列要素数(size 関数で自動的に計算)
  real :: dx  ! x 方向の格子点間隔
  real :: dy  ! y 方向の格子点間隔

  nx=size(x)
  ny=size(y)
  dx=x(2)-x(1)
  dy=y(2)-y(1)

  if(present(undeff))then
     do i=1,ny
        call grad_1d(x, u(:,i), valx(:,i), undeff)
     end do

     do i=1,nx
        call grad_1d(y, u(i,:), valy(i,:), undeff)
     end do
  else
     do i=1,ny
        call grad_1d(x, u(:,i), valx(:,i) )
     end do

     do i=1,nx
        call grad_1d(y, u(i,:), valy(i,:) )
     end do
  end if

end subroutine grad_vec_2d
Subroutine :
x(:) :real, intent(in)
: x 方向の座標変数 [m]
y(:) :real, intent(in)
: y 方向の座標変数 [m]
z(:) :real, intent(in)
: z 方向の座標変数 [m]
u(size(x),size(y),size(z)) :real, intent(in)
: 勾配をとる 2 次元スカラー成分
valx(size(x),size(y),size(z)) :real, intent(inout)
: 計算された y 方向の 2 次元勾配ベクトル
valy(size(x),size(y),size(z)) :real, intent(inout)
: 計算された y 方向の 2 次元勾配ベクトル
valz(size(x),size(y),size(z)) :real, intent(inout)
: 計算された z 方向の 2 次元勾配ベクトル
undeff :real, intent(in), optional
hx(size(x),size(y),size(z)) :real, intent(in), optional
: x 方向のスケール因子
hy(size(x),size(y),size(z)) :real, intent(in), optional
: y 方向のスケール因子
hz(size(x),size(y),size(z)) :real, intent(in), optional
: z 方向のスケール因子

$nabla p =left(frac{partial p}{partial x} ,\; frac{partial p}{partial y} ,\; frac{partial p}{partial z} right) $ を 2 次の中央差分近似で書き換えると, 点 $(i,j,k)$ での勾配は $left(frac{p_{i+1,j,k}-p_{i-1,j,k}}{2dx} ,\; frac{p_{i,j+1,k}-p_{i,j-1,k}}{2dy} ,\; frac{p_{i,j,k+1}-p_{i,j,k-1}}{2dz} right) $ とできる. これを用いて 3 次元勾配を計算. データ点が足りない隅の領域では, 1 次の差分近似で計算するので, 少し精度が 落ちる. 基本的に等間隔格子のデータのみしか計算できないので, 不等間隔のデータは 等間隔に内挿なりを行ってからこのルーチンを呼び出す必要あり. 先に用いた 1 次元勾配計算ルーチンを 3 回呼び出すことにしている.


[Source]

subroutine grad_vec_3d( x, y, z, u, valx, valy, valz, undeff, hx, hy, hz )
  !-------------------------------------------------------------------
  ! $\nabla p =\left(\frac{\partial p}{\partial x} ,\; \frac{\partial p}{\partial y} ,\; \frac{\partial p}{\partial z} \right) $ を
  ! 2 次の中央差分近似で書き換えると, 点 $(i,j,k)$ での勾配は
  ! $\left(\frac{p_{i+1,j,k}-p_{i-1,j,k}}{2dx} ,\; \frac{p_{i,j+1,k}-p_{i,j-1,k}}{2dy} ,\; \frac{p_{i,j,k+1}-p_{i,j,k-1}}{2dz} \right) $
  ! とできる. これを用いて 3 次元勾配を計算.
  ! データ点が足りない隅の領域では, 1 次の差分近似で計算するので, 少し精度が
  ! 落ちる.
  ! 基本的に等間隔格子のデータのみしか計算できないので, 不等間隔のデータは
  ! 等間隔に内挿なりを行ってからこのルーチンを呼び出す必要あり.
  ! 先に用いた 1 次元勾配計算ルーチンを 3 回呼び出すことにしている.
  !-------------------------------------------------------------------
  implicit none
  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))  ! 勾配をとる 2 次元スカラー成分
  real, intent(inout) :: valx(size(x),size(y),size(z))  ! 計算された y 方向の 2 次元勾配ベクトル
  real, intent(inout) :: valy(size(x),size(y),size(z))  ! 計算された y 方向の 2 次元勾配ベクトル
  real, intent(inout) :: valz(size(x),size(y),size(z))  ! 計算された z 方向の 2 次元勾配ベクトル
  real, intent(in), optional :: undeff
  real, intent(in), optional :: hx(size(x),size(y),size(z))  ! x 方向のスケール因子
  real, intent(in), optional :: hy(size(x),size(y),size(z))  ! y 方向のスケール因子
  real, intent(in), optional :: hz(size(x),size(y),size(z))  ! z 方向のスケール因子
  integer :: i   ! イタレーション用添字
  integer :: j   ! イタレーション用添字
  integer :: k   ! イタレーション用添字
  integer :: nx  ! x の配列要素数(size 関数で自動的に計算)
  integer :: ny  ! y の配列要素数(size 関数で自動的に計算)
  integer :: nz  ! z の配列要素数(size 関数で自動的に計算)
  real :: dx  ! x 方向の格子点間隔
  real :: dy  ! y 方向の格子点間隔
  real :: dz  ! z 方向の格子点間隔

  nx=size(x)
  ny=size(y)
  nz=size(z)
  dx=x(2)-x(1)
  dy=y(2)-y(1)
  dz=z(2)-z(1)

  if(present(undeff))then
     do i=1,nz
        do j=1,ny
           call grad_1d(x, u(:,j,i), valx(:,j,i),undeff)
        end do
     end do

     do i=1,nz
        do j=1,nx
           call grad_1d(y, u(j,:,i), valy(j,:,i),undeff)
        end do
     end do

     do i=1,ny
        do j=1,nx
           call grad_1d(z, u(j,i,:), valz(j,i,:),undeff)
        end do
     end do
  else
     do i=1,nz
        do j=1,ny
           call grad_1d(x, u(:,j,i), valx(:,j,i))
        end do
     end do

     do i=1,nz
        do j=1,nx
           call grad_1d(y, u(j,:,i), valy(j,:,i))
        end do
     end do

     do i=1,ny
        do j=1,nx
           call grad_1d(z, u(j,i,:), valz(j,i,:))
        end do
     end do
  end if

end subroutine grad_vec_3d
Subroutine :
x(:) :real, intent(in)
: x 方向の座標変数 [m]
u(size(x)) :real, intent(in)
: 上の空間配列に対応する 1 次元スカラー値
val(size(x)) :real, intent(inout)
: スカラー値の x 方向の勾配
undef :real, intent(in), optional
hx(size(x)) :real, intent(in), optional
: x 方向のスケール因子

1 次元のスカラー変数のラプラシアンを計算する $frac{partial ^2p}{partial x^2} $ を 2 次の中央差分近似で書き換えると, 点 $(i)$ での勾配は $frac{p_{i+1}+p_{i-1}-2p_i}{dx^2} $ とできる. これを用いて 1 次元ラプラシアンを計算. データ点が足りない隅の領域では, undef を定義する. option で undef が定義されていない場合は, 0.0 を代入する. 基本的に等間隔格子のデータのみしか計算できないので, 不等間隔のデータは 等間隔に内挿なりを行ってからこのルーチンを呼び出す必要あり.

[Source]

subroutine laplacian_1d( x, u, val, undef, hx )
! 1 次元のスカラー変数のラプラシアンを計算する
! $\frac{\partial ^2p}{\partial x^2} $ を
! 2 次の中央差分近似で書き換えると, 点 $(i)$ での勾配は
! $\frac{p_{i+1}+p_{i-1}-2p_i}{dx^2} $
! とできる. これを用いて 1 次元ラプラシアンを計算.
! データ点が足りない隅の領域では, undef を定義する.
! option で undef が定義されていない場合は, 0.0 を代入する.
! 基本的に等間隔格子のデータのみしか計算できないので, 不等間隔のデータは
! 等間隔に内挿なりを行ってからこのルーチンを呼び出す必要あり.
  implicit none
  real, intent(in) :: x(:)  ! x 方向の座標変数 [m]
  real, intent(in) :: u(size(x))  ! 上の空間配列に対応する 1 次元スカラー値
  real, intent(inout) :: val(size(x))  ! スカラー値の x 方向の勾配
  real, intent(in), optional :: undef
  real, intent(in), optional :: hx(size(x))  ! x 方向のスケール因子
  integer :: i  ! イタレーション用添字
  integer :: nx  ! 配列要素数
  real :: dx  ! 格子間隔 [m]

  nx=size(x)
  dx=x(2)-x(1)

  if(present(undef))then
     do i=2,nx-1
        if(u(i+1)==undef.or.u(i-1)==undef.or.u(i)==undef)then
           val(i)=undef
        else
           val(i)=(u(i+1)+u(i-1)-2.0*u(i))/(dx*dx)
        end if
     end do
!-- データ数のない両端の処理 ---
!     if(u(1)==undef.or.u(2)==undef)then
!        val(1)=undef
!     else
!        val(1)=(u(2)-u(1))/dx
!     end if
!     if(u(nx)==undef.or.u(nx-1)==undef)then
!        val(nx)=undef
!     else
!        val(nx)=(u(nx)-u(nx-1))/dx
!     end if
     val(1)=undef
     val(nx)=undef
  else
     do i=2,nx-1
        val(i)=(u(i+1)+u(i-1)-2.0*u(i))/(dx*dx)
     end do
!-- データ数のない両端の処理 ---
!     val(1)=(u(2)-u(1))/dx
!     val(nx)=(u(nx)-u(nx-1))/dx
     val(1)=0.0
     val(nx)=0.0
  end if

end subroutine laplacian_1d
Subroutine :
xp :real, intent(in)
: 中心位置座標の x 成分
yp :real, intent(in)
: 中心位置座標の y 成分
zp :real, intent(in)
: 中心位置座標の z 成分
x(:) :real, intent(in)
: x 方向の位置座標
y(:) :real, intent(in)
: y 方向の位置座標
z(:) :real, intent(in)
: z 方向の位置座標
rad(size(x),size(y),size(z)) :real, intent(inout)
: 距離

ある位置からの距離を計算するルーチン 配列数を調整することにより, 2 次元での計算も可能.

[Source]

subroutine radius(xp,yp,zp,x,y,z,rad)
  ! ある位置からの距離を計算するルーチン
  ! 配列数を調整することにより, 2 次元での計算も可能.
  implicit none
  real, intent(in) :: xp  ! 中心位置座標の x 成分
  real, intent(in) :: yp  ! 中心位置座標の y 成分
  real, intent(in) :: zp  ! 中心位置座標の z 成分
  real, intent(in) :: x(:)  ! x 方向の位置座標
  real, intent(in) :: y(:)  ! y 方向の位置座標
  real, intent(in) :: z(:)  ! z 方向の位置座標
  real, intent(inout) :: rad(size(x),size(y),size(z))  ! 距離
  integer :: i, j, k, nx, ny, nz

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

!$omp parallel do shared(rad,x,y,z,xp,yp,zp) private(i,j,k)
  do k=1,nz
     do j=1,ny
        do i=1,nx
           rad(i,j,k)=sqrt((x(i)-xp)**2+(y(j)-yp)**2+(z(k)-zp)**2)
        end do
     end do
  end do
!$omp end parallel do


end subroutine radius
Subroutine :
x(:) :real, intent(in)
: x 方向の空間座標 [m]
y(:) :real, intent(in)
: y 方向の空間座標 [m]
z(:) :real, intent(in)
: z 方向の空間座標 [m]
u(size(x),size(y),size(z)) :real, intent(in)
: x に対応する方向の 2 次元ベクトル成分
v(size(x),size(y),size(z)) :real, intent(in)
: y に対応する方向の 2 次元ベクトル成分
w(size(x),size(y),size(z)) :real, intent(in)
: y に対応する方向の 2 次元ベクトル成分
zeta(size(x),size(y),size(z)) :real, intent(inout)
: 渦度ベクトル x 成分
eta(size(x),size(y),size(z)) :real, intent(inout)
: 渦度ベクトル y 成分
xi(size(x),size(y),size(z)) :real, intent(inout)
: 渦度ベクトル z 成分
undeff :real, intent(in), optional
hx(size(x),size(y),size(z)) :real, intent(in), optional
: x 方向のスケール因子
hy(size(x),size(y),size(z)) :real, intent(in), optional
: y 方向のスケール因子
hz(size(x),size(y),size(z)) :real, intent(in), optional
: z 方向のスケール因子

3 次元渦度を計算する. 引数の順番は右手系で x, y, z のデカルト座標系, それに対応するベクトル成分 u, v, w を代入すると, それに対応した渦度ベクトル成分 zeta, eta, xi が計算される. curl の使いまわし. div_3d でも div を使い回したが, rotate は各面に垂直なベクトルとして 渦度ベクトルを計算するので, 境界で面と境界線を分離して計算する 必要はない. その面に垂直な渦度ベクトルは curl で境界まで考慮して計算可能である. 将来, 任意の格子座標系に対応するため, 後ろにスケーリングファクターオプションを 追加したが, 現在は使用していない.

[Source]

subroutine rotate( x, y, z, u, v, w, zeta, eta, xi, undeff, hx, hy, hz )
! 3 次元渦度を計算する.
! 引数の順番は右手系で x, y, z のデカルト座標系,
! それに対応するベクトル成分 u, v, w を代入すると,
! それに対応した渦度ベクトル成分 zeta, eta, xi が計算される.
! curl の使いまわし.
! div_3d でも div を使い回したが, rotate は各面に垂直なベクトルとして
! 渦度ベクトルを計算するので, 境界で面と境界線を分離して計算する
! 必要はない.
! その面に垂直な渦度ベクトルは curl で境界まで考慮して計算可能である.
! 将来, 任意の格子座標系に対応するため, 後ろにスケーリングファクターオプションを
! 追加したが, 現在は使用していない.
  implicit none
  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))  ! x に対応する方向の 2 次元ベクトル成分
  real, intent(in) :: v(size(x),size(y),size(z))  ! y に対応する方向の 2 次元ベクトル成分
  real, intent(in) :: w(size(x),size(y),size(z))  ! y に対応する方向の 2 次元ベクトル成分
  real, intent(inout) :: zeta(size(x),size(y),size(z))  ! 渦度ベクトル x 成分
  real, intent(inout) :: eta(size(x),size(y),size(z))  ! 渦度ベクトル y 成分
  real, intent(inout) :: xi(size(x),size(y),size(z))  ! 渦度ベクトル z 成分
  real, intent(in), optional :: undeff
  real, intent(in), optional :: hx(size(x),size(y),size(z))  ! x 方向のスケール因子
  real, intent(in), optional :: hy(size(x),size(y),size(z))  ! y 方向のスケール因子
  real, intent(in), optional :: hz(size(x),size(y),size(z))  ! z 方向のスケール因子
  integer :: i   ! イタレーション用添字
  integer :: j   ! イタレーション用添字
  integer :: k   ! イタレーション用添字
  integer :: nx  ! 空間配列要素数 1 次元目
  integer :: ny  ! 空間配列要素数 2 次元目
  integer :: nz  ! 空間配列要素数 3 次元目
  real :: dx  ! 1 次元目を微分する格子間隔 [m]
  real :: dy  ! 2 次元目を微分する格子間隔 [m]
  real :: dz  ! 3 次元目を微分する格子間隔 [m]

  dx=x(2)-x(1)
  dy=y(2)-y(1)
  dz=z(2)-z(1)
  nx=size(x)
  ny=size(y)
  nz=size(z)

  if(present(undeff))then
!$omp parallel default(shared)
!$omp do private(i)
! x 軸に垂直な面上での x 方向の渦度ベクトルを 3 次元全域で計算.
     do i=1,nx
        call curl( y, z, v(i,:,:), w(i,:,:), zeta(i,:,:), undeff )
     end do

!$omp end do
!$omp end parallel
  else
!$omp parallel default(shared)
!$omp do private(i)
! x 軸に垂直な面上での x 方向の渦度ベクトルを 3 次元全域で計算.
     do i=1,nx
        call curl( y, z, v(i,:,:), w(i,:,:), zeta(i,:,:) )
     end do
!$omp end do
!$omp end parallel

  end if

! y 軸に垂直な面上での y 方向の渦度ベクトルを 3 次元全域で計算.

  if(present(undeff))then
!$omp parallel default(shared)
!$omp do private(j)
     do j=1,ny
        call curl( x, z, u(:,j,:), w(:,j,:), eta(:,j,:), undeff, ord=.false. )
     end do
!$omp end do
!$omp end parallel

  else
!$omp parallel default(shared)
!$omp do private(j)
     do j=1,ny
        call curl( x, z, u(:,j,:), w(:,j,:), eta(:,j,:), ord=.false. )
     end do
!$omp end do
!$omp end parallel

  end if

! z 軸に垂直な面上での z 方向の渦度ベクトルを 3 次元全域で計算.

  if(present(undeff))then
!$omp parallel default(shared)
!$omp do private(k)
     do k=1,nz
        call curl( x, y, u(:,:,k), v(:,:,k), xi(:,:,k), undeff )
     end do

!$omp end do
!$omp end parallel
  else
!$omp parallel default(shared)
!$omp do private(k)
     do k=1,nz
        call curl( x, y, u(:,:,k), v(:,:,k), xi(:,:,k) )
     end do

!$omp end do
!$omp end parallel
  end if

end subroutine
Subroutine :
x(:,:,:) :real, intent(in)
: x 方向のベクトル成分
y(size(x,1),size(x,2),size(x,3)) :real, intent(in)
: y 方向のベクトル成分
z(size(x,1),size(x,2),size(x,3)) :real, intent(in)
: z 方向のベクトル成分
u(size(x,1),size(x,2),size(x,3)) :real, intent(in)
: x 方向のベクトル成分
v(size(x,1),size(x,2),size(x,3)) :real, intent(in)
: y 方向のベクトル成分
w(size(x,1),size(x,2),size(x,3)) :real, intent(in)
: z 方向のベクトル成分
vecx(size(x,1),size(x,2),size(x,3)) :real, intent(inout)
: 外積の x 成分
vecy(size(x,1),size(x,2),size(x,3)) :real, intent(inout)
: 外積の y 成分
vecz(size(x,1),size(x,2),size(x,3)) :real, intent(inout)
: 外積の z 成分
undeff :real, intent(in), optional

2ベクトルの外積計算ルーチン 配列の要素数を工夫することで 2 次元外積も計算可能

[Source]

subroutine vec_prod(x,y,z,u,v,w,vecx,vecy,vecz,undeff)
  ! 2ベクトルの外積計算ルーチン
  ! 配列の要素数を工夫することで 2 次元外積も計算可能
  implicit none
  real, intent(in) :: x(:,:,:)  ! x 方向のベクトル成分
  real, intent(in) :: y(size(x,1),size(x,2),size(x,3))  ! y 方向のベクトル成分
  real, intent(in) :: z(size(x,1),size(x,2),size(x,3))  ! z 方向のベクトル成分
  real, intent(in) :: u(size(x,1),size(x,2),size(x,3))  ! x 方向のベクトル成分
  real, intent(in) :: v(size(x,1),size(x,2),size(x,3))  ! y 方向のベクトル成分
  real, intent(in) :: w(size(x,1),size(x,2),size(x,3))  ! z 方向のベクトル成分
  real, intent(inout) :: vecx(size(x,1),size(x,2),size(x,3))  ! 外積の x 成分
  real, intent(inout) :: vecy(size(x,1),size(x,2),size(x,3))  ! 外積の y 成分
  real, intent(inout) :: vecz(size(x,1),size(x,2),size(x,3))  ! 外積の z 成分
  real, intent(in), optional :: undeff
  integer :: i, j, k, nx, ny, nz

  nx=size(x,1)
  ny=size(x,2)
  nz=size(x,3)

  if(present(undeff))then
!$omp parallel do shared(vecx,vecy,vecz,x,y,z,u,v,w) private(i,j,k)
     do k=1,nz
        do j=1,ny
           do i=1,nx
              if(x(i,j,k)==undeff.or.u(i,j,k)==undeff.or.y(i,j,k)==undeff.or. v(i,j,k)==undeff.or.z(i,j,k)==undeff.or.w(i,j,k)==undeff)then
                 vecx(i,j,k)=undeff
                 vecy(i,j,k)=undeff
                 vecz(i,j,k)=undeff
              else
                 vecx(i,j,k)=y(i,j,k)*w(i,j,k)-z(i,j,k)*v(i,j,k)
                 vecy(i,j,k)=z(i,j,k)*u(i,j,k)-x(i,j,k)*w(i,j,k)
                 vecz(i,j,k)=x(i,j,k)*v(i,j,k)-y(i,j,k)*u(i,j,k)
              end if
           end do
        end do
     end do
!$omp end parallel do

  else

!$omp parallel do shared(vecx,vecy,vecz,x,y,z,u,v,w) private(i,j,k)
     do k=1,nz
        do j=1,ny
           do i=1,nx
              vecx(i,j,k)=y(i,j,k)*w(i,j,k)-z(i,j,k)*v(i,j,k)
              vecy(i,j,k)=z(i,j,k)*u(i,j,k)-x(i,j,k)*w(i,j,k)
              vecz(i,j,k)=x(i,j,k)*v(i,j,k)-y(i,j,k)*u(i,j,k)
           end do
        end do
     end do
!$omp end parallel do
  end if

end subroutine vec_prod