Class analy
In: analy.f90

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

Methods

Public Instance methods

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

任意の点まわりの絶対角運動量を計算するルーチン 主目的は台風の中心軸を中心に鉛直軸まわりの角運動量を計算することを目的とする. $M=rv+dfrac{fr^2}{2} ,quad r=中心軸からの距離, quad v=風速の同位角成分$ 位置と風速に緯度の変換を与えれば, 全球での自転軸まわりの角運動量を 計算することも可能. ベクトルの外積計算ルーチン vec_prod を用いることで極座標でも計算可能.

[Source]

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

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

  tmp=0.0
  call vec_prod(nx,ny,1,xxp,yyp,tmp,u,v,tmp,tmp,tmp,mome)
  call radius(nx,ny,1,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 :
nx :integer, intent(in)
: x 方向の配列数
ny :integer, intent(in)
: y 方向の配列数
xp :real, intent(in)
: 回転軸の x 位置座標
yp :real, intent(in)
: 回転軸の y 位置座標
x(nx) :real, intent(in)
: x 方向の位置座標
y(ny) :real, intent(in)
: y 方向の位置座標
u(nx,ny) :real, intent(in)
: 位置 i,j での風速の 1 成分
v(nx,ny) :real, intent(in)
: 位置 i,j での風速の 1 成分
mome(nx,ny) :real, intent(inout)
: 回転軸まわりの相対角運動量

任意の点まわりの相対角運動量を計算するルーチン 本当は 3 次元ベクトルで計算されるが, 気象学では 3 次元量はあまり需要がない であろうという判断から, ある回転軸まわりの角運動量成分のみを 計算することにしている. 主目的は台風の中心軸を中心に鉛直軸まわりの角運動量を計算することを目的とする. $M=rv,quad r=中心軸からの距離, quad v=風速の同位角成分$ 位置と風速に緯度の変換を与えれば, 全球での自転軸まわりの角運動量を 計算することも可能. ベクトルの外積計算ルーチン vec_prod を用いることで極座標でも計算可能.

[Source]

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

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

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

end subroutine Rangular_moment
Subroutine :
nx :integer, intent(in)
: x 方向の配列数
ny :integer, intent(in)
: y 方向の配列数
nz :integer, intent(in)
: z 方向の配列数
x(nx,ny,nz) :real, intent(in)
: x 方向のベクトル成分
y(nx,ny,nz) :real, intent(in)
: y 方向のベクトル成分
z(nx,ny,nz) :real, intent(in)
: z 方向のベクトル成分
dis(nx,ny,nz) :real, intent(inout)
: 各点での絶対値ベクトル

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

[Source]

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

!$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)
: 空間配列 1 次元目
y(:) :real, intent(in)
: 空間配列 2 次元目
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(out)
: 2次元渦度値

— x,y と u,v は偶置換の向きに配置すれば, 任意の2次元渦度が計算可能 —


$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)
!-- x,y と u,v は偶置換の向きに配置すれば, 任意の2次元渦度が計算可能 ---
  !-------------------------------------------------------------------
  ! $\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(:)   ! 空間配列 1 次元目
  real, intent(in) :: y(:)   ! 空間配列 2 次元目
  real, intent(in) :: u(size(x),size(y))  ! x に対応する方向の 2 次元ベクトル成分
  real, intent(in) :: v(size(x),size(y))  ! y に対応する方向の 2 次元ベクトル成分
  real, intent(out) :: val(size(x),size(y))  ! 2次元渦度値
  integer :: nx  ! x の配列要素数(size 関数で自動的に計算)
  integer :: ny  ! y の配列要素数(size 関数で自動的に計算)
  integer :: i   ! イタレーション用添字
  integer :: j   ! イタレーション用添字
  real :: dx  ! x 方向の格子点間隔
  real :: dy  ! y 方向の格子点間隔

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

  dx=x(2)-x(1)
  dy=y(2)-y(1)

!$omp parallel do shared(val,u,v,dx,dy) private(i,j)
  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
!$omp end parallel 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 subroutine curl
Subroutine :
x(:) :real, intent(in)
: 空間配列 1 次元目
y(:) :real, intent(in)
: 空間配列 2 次元目
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(out)
: 2次元発散値

2次元発散計算ルーチン — x,y と u,v は偶置換の向きに配置すれば, 任意の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)   ! 2次元発散計算ルーチン
!-- x,y と u,v は偶置換の向きに配置すれば, 任意の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(:)   ! 空間配列 1 次元目
  real, intent(in) :: y(:)   ! 空間配列 2 次元目
  real, intent(in) :: u(size(x),size(y))  ! x に対応する方向の 2 次元ベクトル成分
  real, intent(in) :: v(size(x),size(y))  ! y に対応する方向の 2 次元ベクトル成分
  real, intent(out) :: val(size(x),size(y))  ! 2次元発散値
  integer :: nx  ! x の配列要素数(size 関数で自動的に計算)
  integer :: ny  ! y の配列要素数(size 関数で自動的に計算)
  integer :: i   ! イタレーション用添字
  integer :: j   ! イタレーション用添字
  real :: dx  ! x 方向の格子点間隔
  real :: dy  ! y 方向の格子点間隔

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

  dx=x(2)-x(1)
  dy=y(2)-y(1)

!$omp parallel do shared(val,u,v,dx,dy) private(i,j)
  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
!$omp end parallel 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 subroutine div
Subroutine :
nx :integer, intent(in)
: x 方向の配列数
ny :integer, intent(in)
: y 方向の配列数
nz :integer, intent(in)
: z 方向の配列数
x(nx,ny,nz) :real, intent(in)
: x 方向のベクトル成分
y(nx,ny,nz) :real, intent(in)
: y 方向のベクトル成分
z(nx,ny,nz) :real, intent(in)
: z 方向のベクトル成分
u(nx,ny,nz) :real, intent(in)
: x 方向のベクトル成分
v(nx,ny,nz) :real, intent(in)
: y 方向のベクトル成分
w(nx,ny,nz) :real, intent(in)
: z 方向のベクトル成分
dot(nx,ny,nz) :real, intent(inout)
: 内積

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

[Source]

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

!$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 subroutine dot_prod
Subroutine :
x(:) :real, intent(in)
: 空間配列 1 次元
u(size(x)) :real, intent(in)
: 上の空間配列に対応する 1 次元スカラー値
val(size(x)) :real, intent(inout)
: スカラー値の x 方向の勾配

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


[Source]

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

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

  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 subroutine grad_1d
Subroutine :
x(:) :real, intent(in)
: 空間配列 1 次元目
y(:) :real, intent(in)
: 空間配列 2 次元目
u(size(x),size(y)) :real, intent(in)
: 勾配をとる 2 次元スカラー成分
valx(size(x),size(y)) :real, intent(out)
: 計算された x 方向の 2 次元勾配ベクトル
valy(size(x),size(y)) :real, intent(out)
: 計算された x 方向の 2 次元勾配ベクトル

$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 )
  !-------------------------------------------------------------------
  ! $\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(:)   ! 空間配列 1 次元目
  real, intent(in) :: y(:)   ! 空間配列 2 次元目
  real, intent(in) :: u(size(x),size(y))  ! 勾配をとる 2 次元スカラー成分
  real, intent(out) :: valx(size(x),size(y))  ! 計算された x 方向の 2 次元勾配ベクトル
  real, intent(out) :: valy(size(x),size(y))  ! 計算された x 方向の 2 次元勾配ベクトル
  integer :: nx  ! x の配列要素数(size 関数で自動的に計算)
  integer :: ny  ! y の配列要素数(size 関数で自動的に計算)
  integer :: i   ! イタレーション用添字
  integer :: j   ! イタレーション用添字
  real :: dx  ! x 方向の格子点間隔
  real :: dy  ! y 方向の格子点間隔

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

  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 subroutine grad_vec_2d
Subroutine :
x(:) :real, intent(in)
: 空間配列 1 次元目
y(:) :real, intent(in)
: 空間配列 2 次元目
z(:) :real, intent(in)
: 空間配列 3 次元目
u(size(x),size(y),size(z)) :real, intent(in)
: 勾配をとる 3 次元スカラー成分
valx(size(x),size(y),size(z)) :real, intent(inout)
: 計算された x 方向の 3 次元勾配ベクトル
valy(size(x),size(y),size(z)) :real, intent(inout)
: 計算された y 方向の 3 次元勾配ベクトル
valz(size(x),size(y),size(z)) :real, intent(inout)
: 計算された z 方向の 3 次元勾配ベクトル

$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 )
  !-------------------------------------------------------------------
  ! $\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(:)   ! 空間配列 1 次元目
  real, intent(in) :: y(:)   ! 空間配列 2 次元目
  real, intent(in) :: z(:)   ! 空間配列 3 次元目
  real, intent(in) :: u(size(x),size(y),size(z))  ! 勾配をとる 3 次元スカラー成分
  real, intent(inout) :: valx(size(x),size(y),size(z))  ! 計算された x 方向の 3 次元勾配ベクトル
  real, intent(inout) :: valy(size(x),size(y),size(z))  ! 計算された y 方向の 3 次元勾配ベクトル
  real, intent(inout) :: valz(size(x),size(y),size(z))  ! 計算された z 方向の 3 次元勾配ベクトル
  integer :: nx  ! x の配列要素数(size 関数で自動的に計算)
  integer :: ny  ! y の配列要素数(size 関数で自動的に計算)
  integer :: nz  ! z の配列要素数(size 関数で自動的に計算)
  integer :: i   ! イタレーション用添字
  integer :: j   ! イタレーション用添字
  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)

  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 subroutine grad_vec_3d
Subroutine :
nx :integer, intent(in)
: x 方向の配列数
ny :integer, intent(in)
: y 方向の配列数
nz :integer, intent(in)
: z 方向の配列数
xp :real, intent(in)
: 中心位置座標の x 成分
yp :real, intent(in)
: 中心位置座標の y 成分
zp :real, intent(in)
: 中心位置座標の z 成分
x(nx) :real, intent(in)
: x 方向の位置座標
y(ny) :real, intent(in)
: y 方向の位置座標
z(nz) :real, intent(in)
: z 方向の位置座標
rad(nx,ny,nz) :real, intent(inout)

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

[Source]

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

!$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 :
nx :integer, intent(in)
: x 方向の配列数
ny :integer, intent(in)
: y 方向の配列数
nz :integer, intent(in)
: z 方向の配列数
x(nx,ny,nz) :real, intent(in)
: x 方向のベクトル成分
y(nx,ny,nz) :real, intent(in)
: y 方向のベクトル成分
z(nx,ny,nz) :real, intent(in)
: z 方向のベクトル成分
u(nx,ny,nz) :real, intent(in)
: x 方向のベクトル成分
v(nx,ny,nz) :real, intent(in)
: y 方向のベクトル成分
w(nx,ny,nz) :real, intent(in)
: z 方向のベクトル成分
vecx(nx,ny,nz) :real, intent(inout)
: 外積の x 成分
vecy(nx,ny,nz) :real, intent(inout)
: 外積の y 成分
vecz(nx,ny,nz) :real, intent(inout)
: 外積の z 成分

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

[Source]

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

!$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 subroutine vec_prod