Class Ellip_Slv
In: ellip_slv.f90

代数演算を用いて楕円型偏微分方程式を解くモジュール

Methods

Included Modules

omp_lib

Public Instance methods

Subroutine :
x(:) :real, intent(in)
: 領域の横座標
y(:) :real, intent(in)
: 領域の縦座標
rho(size(x),size(y)) :real, intent(in)
: ポアソン方程式の強制項 rho =0 でラプラス方程式も求積可能
eps :real, intent(in)
: 収束条件
boundary :character(4), intent(in)
: 境界条件 4 文字で各辺の境界条件を与える. 1 文字目 : x 下端, 2 文字目 : y 左端, 3 文字目 : x 上端, 4 文字目 : y 右端 boundary は 1 : 固定端境界, 2 : 自由端境界, 3 : 周期境界
psi(size(x),size(y)) :real, intent(inout)
: ポアソン方程式の解
bound_opt(size(x),size(y)) :real, intent(in), optional
: 境界での強制 ノイマン境界の場合 : フラックス値
a(size(x),size(y)) :real, intent(in), optional
: 各微分項の係数
b(size(x),size(y)) :real, intent(in), optional
: 各微分項の係数
c(size(x),size(y)) :real, intent(in), optional
: 各微分項の係数
d(size(x),size(y)) :real, intent(in), optional
: 各微分項の係数
e(size(x),size(y)) :real, intent(in), optional
: 各微分項の係数
f(size(x),size(y)) :real, intent(in), optional
: 各微分項の係数
undef :real, intent(in), optional
: 未定義値
inner_bound(size(x),size(y)) :integer, intent(in), optional
: 内部領域の境界. 値に応じてその格子点で境界値計算 1 = 固定端境界, 10 = 境界の内側. 2 = y 方向自由端境界 (フラックスは上向き) -2 = y 方向自由端境界 (フラックスは下向き) 4 = x 方向自由端境界 (フラックスは右向き) -4 = x 方向自由端境界 (フラックスは左向き) 3 = 周期境界 8 = |_, ~| で両方とも自由境界条件 -8 = |~, _| で両方とも自由境界条件 この引数が与えられなければ全領域を計算する. 境界の内側格子点 (10) は反復計算を行わず, undef で設定された値もしくはゼロが入る. このときの境界値は bound_opt の値が用いられる.
init_flag :logical, intent(in), optional
: psi の値をゼロで初期化するか. .true. = 初期化する. .false. = 初期化しない. デフォルトでは初期化する.
accel :real, intent(in), optional
: SOR の加速係数 (0 < accel < 2) デフォルト = 1
ln :integer, intent(in), optional
: 反復回数 この値が与えられるとき, eps の値に無関係に ln 回ループさせる.

ガウス=ザイデル法による楕円型方程式の求積 各オプション配列は, ポアソン系の各微分項の値. デフォルトはゼロで設定される. $$a\dfrac{partial ^2\psi}{partial x^2} +b\dfrac{partial ^2\psi}{partial x\partial y} +c\dfrac{partial ^2\psi}{partial y^2} +d\dfrac{partial psi}{partial x} +e\dfrac{partial psi}{partial y} +f\psi =rho $$ の各係数に対応している.

[Source]

subroutine Ellip_GauSei_2d( x, y, rho, eps, boundary, psi, bound_opt, a, b, c, d, e, f, undef, inner_bound, init_flag, accel, ln )
! ガウス=ザイデル法による楕円型方程式の求積
! 各オプション配列は, ポアソン系の各微分項の値. デフォルトはゼロで設定される.
! $$a\dfrac{\partial ^2\psi}{\partial x^2} +b\dfrac{\partial ^2\psi}{\partial x\partial y} +c\dfrac{\partial ^2\psi}{\partial y^2} +d\dfrac{\partial \psi}{\partial x} +e\dfrac{\partial \psi}{\partial y} +f\psi =\rho $$
! の各係数に対応している.
  implicit none
  real, intent(in) :: x(:)  ! 領域の横座標
  real, intent(in) :: y(:)  ! 領域の縦座標
  real, intent(in) :: rho(size(x),size(y))  ! ポアソン方程式の強制項
                   ! rho =0 でラプラス方程式も求積可能
  real, intent(in) :: eps  ! 収束条件
  character(4), intent(in) :: boundary  ! 境界条件
                ! 4 文字で各辺の境界条件を与える.
                ! 1 文字目 : x 下端, 2 文字目 : y 左端, 3 文字目 : x 上端,
                ! 4 文字目 : y 右端
                ! boundary は 1 : 固定端境界, 2 : 自由端境界, 3 : 周期境界
  real, intent(in), optional :: bound_opt(size(x),size(y))  ! 境界での強制
                             ! ノイマン境界の場合 : フラックス値
  real, intent(in), optional :: a(size(x),size(y))  ! 各微分項の係数
  real, intent(in), optional :: b(size(x),size(y))  ! 各微分項の係数
  real, intent(in), optional :: c(size(x),size(y))  ! 各微分項の係数
  real, intent(in), optional :: d(size(x),size(y))  ! 各微分項の係数
  real, intent(in), optional :: e(size(x),size(y))  ! 各微分項の係数
  real, intent(in), optional :: f(size(x),size(y))  ! 各微分項の係数
  real, intent(inout) :: psi(size(x),size(y))  ! ポアソン方程式の解
  real, intent(in), optional :: undef  ! 未定義値
  integer, intent(in), optional :: inner_bound(size(x),size(y))
                             ! 内部領域の境界. 値に応じてその格子点で境界値計算
                             ! 1 = 固定端境界, 10 = 境界の内側.
                             ! 2 = y 方向自由端境界 (フラックスは上向き)
                             ! -2 = y 方向自由端境界 (フラックスは下向き)
                             ! 4 = x 方向自由端境界 (フラックスは右向き)
                             ! -4 = x 方向自由端境界 (フラックスは左向き)
                             ! 3 = 周期境界
                             ! 8 = |_, ~| で両方とも自由境界条件
                             ! -8 = |~, _| で両方とも自由境界条件
                             ! この引数が与えられなければ全領域を計算する.
                             ! 境界の内側格子点 (10) は反復計算を行わず,
                             ! undef で設定された値もしくはゼロが入る.
                             ! このときの境界値は bound_opt の値が用いられる.
  logical, intent(in), optional :: init_flag  ! psi の値をゼロで初期化するか.
                             ! .true. = 初期化する. .false. = 初期化しない.
                             ! デフォルトでは初期化する.
  real, intent(in), optional :: accel  ! SOR の加速係数 (0 < accel < 2)
                             ! デフォルト = 1
  integer, intent(in), optional :: ln  ! 反復回数
                             ! この値が与えられるとき, eps の値に無関係に
                             ! ln 回ループさせる.
  integer :: i, j, ix, jy, nl, counter
  integer :: nx  ! x 方向の配列要素
  integer :: ny  ! y 方向の配列要素
  integer :: signb  ! 各係数を計算するかどうか
  integer, dimension(size(x),size(y)) :: ib

  real :: defun
  real :: tmp, err, err_max
  real :: tmp_b, accc
  real :: bnd(size(x),size(y))
  real :: dx(size(x)), dy(size(y)), dx2(size(x)), dy2(size(y))
  real, dimension(size(x),size(y)) :: dxdy
  real, dimension(size(x),size(y)) :: at, bt, ct, dt, et, ft
  real, dimension(size(x),size(y)) :: adp, adm, cep, cem, ac, divi

  character(4) :: bound
  logical :: sor_flag
  logical, dimension(size(x),size(y)) :: inner_flag

  bound(1:4)=boundary(1:4)

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

!-- 応答関数の初期化

  if(present(init_flag))then
     if(init_flag.eqv..true.)then
        psi = 0.0
     end if
  else
     psi = 0.0
  end if

!-- 内部境界の判別フラグの設定

  if(present(inner_bound))then
     call set_bound( bound, ib, inner_flag, inner_bound )
  else
     call set_bound( bound, ib, inner_flag )
  end if

!-- 領域・内部境界における境界値の設定

  if(present(bound_opt))then
     call setval_bound( ib, bnd, psi, bound_opt )
  else
     call setval_bound( ib, bnd, psi )
  end if

!-- 未定義値の設定

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

!-- 係数の代入
!-- a, c については, 値が入れられていなければ, 全配列に 1 を代入する.
  if(present(a))then
     call set_coe( at, ext=a )
  else
     call set_coe( at, def=1.0 )
  end if

  if(present(c))then
     call set_coe( ct, ext=c )
  else
     call set_coe( ct, def=1.0 )
  end if

  if(present(b))then
     call set_coe( bt, ext=b )
     signb=1
  else
     call set_coe( bt, def=0.0 )
     signb=0
  end if

  if(present(d))then
     call set_coe( dt, ext=d )
  else
     call set_coe( dt, def=0.0 )
  end if

  if(present(e))then
     call set_coe( et, ext=e )
  else
     call set_coe( et, def=0.0 )
  end if

  if(present(f))then
     call set_coe( ft, ext=f )
  else
     call set_coe( ft, def=0.0 )
  end if

!-- 最高階数における係数チェック. (係数がゼロでないか調べる.)

  call check_coe( at, 0.0 )
  call check_coe( ct, 0.0 )

!-- 加速係数の判定

  if(present(accel))then
     accc=accel
     sor_flag=.true.
  else
     sor_flag=.false.
  end if

!-- ln の処理
  if(present(ln))then
     nl=ln
  else
     nl=0
  end if

!-- 以下で先に格子間隔等の 1 回計算でよいものを求めておく.
!-- これらは 1 方向のみで変化すればよい.
!-- 格子点間隔の計算
  do i=2,nx-1
     dx(i)=(x(i+1)-x(i-1))*0.5
     dx2(i)=dx(i)**2
  end do
  do j=2,ny-1
     dy(j)=(y(j+1)-y(j-1))*0.5
     dy2(j)=dy(j)**2
  end do

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

  do j=1,ny
     do i=1,nx
        dxdy(i,j)=dx(i)*dy(j)
     end do
  end do

!-- ポアソン係数の計算
!-- 実際にポアソン係数が使用されるのは, 境界を除いた 1 回り内側なので,
!-- 計算量削減のため, ループをそのようにしておく.

  ac=0.0
  adp=0.0
  adm=0.0
  cep=0.0
  cem=0.0
  bt=0.0

!-- 最高次数係数 ac の計算
!$omp parallel default(shared)
!$omp do schedule(dynamic) private(i,j)

  do j=2,ny-1
     do i=2,nx-1
        ac(i,j)=1.0/(2.0*(at(i,j)/dx2(i)+ct(i,j)/dy2(j))-ft(i,j))
        adp(i,j)=(at(i,j)/(dx2(i))+0.5*dt(i,j)/dx(i))*ac(i,j)
        adm(i,j)=(at(i,j)/(dx2(i))-0.5*dt(i,j)/dx(i))*ac(i,j)
        cep(i,j)=(ct(i,j)/(dy2(j))+0.5*et(i,j)/dy(j))*ac(i,j)
        cem(i,j)=(ct(i,j)/(dy2(j))-0.5*et(i,j)/dy(j))*ac(i,j)
        bt(i,j)=0.25*bt(i,j)/(dxdy(i,j))*ac(i,j)
     end do
  end do

!$omp end do
!$omp end parallel

  err_max=eps  ! while に入るための便宜的措置
  counter=0

!-- 実際のソルバ ---
  do while(err_max>=eps)
     err_max=0.0

     do j=1,ny
        do i=1,nx

!-- 以降, 反復計算に必要な点の周囲 8 点についてそれぞれ
!-- inner_flag のチェックと同時に適切な値をそれぞれ逐次計算する,
!-- 下で, 各方向の格子計算を行っているが, ib=8,-8 は 4 隅格子のときにしか
!-- case select しない. なぜなら, 真横や上下に隅格子がくることはありえない.
           if(inner_flag(i,j).eqv..false.)then  ! .false. なら領域計算開始
              tmp=-rho(i,j)*ac(i,j)
              tmp=tmp+adp(i,j)*psi(i+1,j) +adm(i,j)*psi(i-1,j) +cep(i,j)*psi(i,j+1) +cem(i,j)*psi(i,j-1)

              if(signb==0)then  ! そもそも bt = 0 なら計算しない.
                 tmp_b=0.0
              else
                 tmp_b=bt(i,j)*(psi(i+1,j+1)+psi(i-1,j-1) -psi(i+1,j-1)-psi(i-1,j+1))
              end if

              tmp=tmp+tmp_b

           else  ! .true. なら境界計算に移行.

              select case (ib(i,j))
              case (1)
                 tmp=bnd(i,j)

              case (2)  ! x 方向にフラックス一定, 上側が参照値
                 tmp=psi(i+1,j)-bnd(i,j)*dx(i)

              case (-2)  ! x 方向にフラックス一定, 下側が参照値
                 tmp=psi(i-1,j)+bnd(i,j)*dx(i)

              case (4)  ! y 方向にフラックス一定, 右側が参照値
                 tmp=psi(i,j+1)-bnd(i,j)*dy(j)

              case (-4)  ! y 方向にフラックス一定, 左側が参照値
                 tmp=psi(i,j-1)+bnd(i,j)*dy(j)

              case (3)  ! 周期境界
                 if(i==1)then
                    ix=nx-1
                 else if(i==nx)then
                    ix=2
                 else
                    ix=i
                 end if
                 if(j==1)then
                    jy=ny-1
                 else if(j==ny)then
                    jy=2
                 else
                    jy=j
                 end if

                 tmp=psi(ix,jy)

              case (7)  ! 両方フラックス一定で内部境界限定.
                 if((ib(i+1,j+1)==10).and.(ib(i+1,j)/=10).and. (ib(i,j+1)/=10))then
                    tmp=0.5*(psi(i-1,j+1)+psi(i+1,j-1) +bnd(i,j+1)*dx(i)+bnd(i+1,j)*dy(j))

                 else if((ib(i-1,j-1)==10).and.(ib(i-1,j)/=10).and. (ib(i,j-1)/=10))then
                    tmp=0.5*(psi(i+1,j-1)+psi(i-1,j+1) -bnd(i,j-1)*dx(i)-bnd(i-1,j)*dy(j))

                 end if

              case (-7)  ! 両方フラックス一定で内部境界限定.
                 if((ib(i-1,j+1)==10).and.(ib(i-1,j)/=10).and. (ib(i,j+1)/=10))then
                    tmp=0.5*(psi(i+1,j+1)+psi(i-1,j-1) -bnd(i,j+1)*dx(i)+bnd(i-1,j)*dy(j))

                 else if((ib(i+1,j-1)==10).and.(ib(i+1,j)/=10).and. (ib(i,j-1)/=10))then
                    tmp=0.5*(psi(i-1,j-1)+psi(i+1,j+1) +bnd(i,j-1)*dx(i)-bnd(i+1,j)*dy(j))

                 end if

              case (8)  ! 両方フラックス一定で左下角か右上角, もしくは内部境界.
                 if(i==1.and.j==1)then  ! -- 評価 1
                    tmp=psi(i+1,j+1)-0.5*(bnd(i+1,j)*dy(j)+bnd(i,j+1)*dx(i))

                 else if(i==nx.and.j==ny)then  ! -- 評価 2
                    tmp=psi(i-1,j-1)+0.5*(bnd(i-1,j)*dy(j)+bnd(i,j-1)*dx(i))

                 else if(ib(i-1,j)==10.and.ib(i,j-1)==10)then
                    ! -- 評価 1 と同じ
                    tmp=psi(i+1,j+1)-0.5*(bnd(i+1,j)*dy(j)+bnd(i,j+1)*dx(i))

                 else if(ib(i+1,j)==10.and.ib(i,j+1)==10)then
                    ! -- 評価 2 と同じ
                    tmp=psi(i-1,j-1)+0.5*(bnd(i-1,j)*dy(j)+bnd(i,j-1)*dx(i))

                 end if

              case (-8)  ! 両方フラックス一定で右下角か左上角
                 if(i==1.and.j==ny)then  ! -- 評価 1
                    tmp=psi(i+1,j-1)+0.5*(bnd(i+1,j)*dy(j)-bnd(i,j-1)*dx(i))

                 else if(i==nx.and.j==1)then  ! -- 評価 2
                    tmp=psi(i-1,j+1)+0.5*(-bnd(i-1,j)*dy(j)+bnd(i,j+1)*dx(i))

                 else if(ib(i-1,j)==10.and.ib(i,j+1)==10)then
                    ! -- 評価 1 と同じ
                    tmp=psi(i+1,j-1)+0.5*(bnd(i+1,j)*dy(j)-bnd(i,j-1)*dx(i))

                 else if(ib(i+1,j)==10.and.ib(i,j-1)==10)then
                    ! -- 評価 2 と同じ
                    tmp=psi(i-1,j+1)+0.5*(-bnd(i-1,j)*dy(j)+bnd(i,j+1)*dx(i))

                 end if
              end select

           end if

           if(sor_flag.eqv..true.)then
              tmp=(1.0-accc)*psi(i,j)+tmp*accc
           end if

           err=abs(tmp-psi(i,j))

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

           psi(i,j)=tmp

        end do
     end do

     if(nl/=0)then
        err_max=eps
        counter=counter+1
        if(counter==nl)then
           exit
        end if
     end if

  end do

!-- 境界の設定

  call calculate_bound( ib, dx, dy, bnd, psi )

!-- 未定義領域には undef を代入する.

  do j=1,ny
     do i=1,nx
        if(ib(i,j)==10)then
           psi(i,j)=defun
        end if
     end do
  end do

end subroutine Ellip_GauSei_2d
Subroutine :
x(:) :real, intent(in)
: 領域の x 座標
y(:) :real, intent(in)
: 領域の y 座標
z(:) :real, intent(in)
: 領域の z 座標
rho(size(x),size(y),size(z)) :real, intent(in)
: ポアソン方程式の強制項 rho =0 でラプラス方程式も求積可能
eps :real, intent(in)
: 収束条件
boundary :character(6), intent(in)
: 境界条件 4 文字で各辺の境界条件を与える. 1 文字目 : x 下端, 2 文字目 : y 左端, 3 文字目 : x 上端, 4 文字目 : y 右端 boundary は 1 : 固定端境界, 2 : 自由端境界, 3 : 周期境界
psi(size(x),size(y),size(z)) :real, intent(inout)
: ポアソン方程式の解
bound_opt(size(x),size(y),size(z)) :real, intent(in), optional
: 境界での強制 ノイマン境界の場合 : フラックス値
xa(size(x),size(y),size(z)) :real, intent(in), optional
: 各微分項の係数
ya(size(x),size(y),size(z)) :real, intent(in), optional
: 各微分項の係数
za(size(x),size(y),size(z)) :real, intent(in), optional
: 各微分項の係数
a(size(x),size(y),size(z)) :real, intent(in), optional
: 各微分項の係数
b(size(x),size(y),size(z)) :real, intent(in), optional
: 各微分項の係数
c(size(x),size(y),size(z)) :real, intent(in), optional
: 各微分項の係数
d(size(x),size(y),size(z)) :real, intent(in), optional
: 各微分項の係数
e(size(x),size(y),size(z)) :real, intent(in), optional
: 各微分項の係数
f(size(x),size(y),size(z)) :real, intent(in), optional
: 各微分項の係数
g(size(x),size(y),size(z)) :real, intent(in), optional
: 各微分項の係数
undef :real, intent(in), optional
: 未定義値
inner_bound(size(x),size(y),size(z)) :integer, intent(in), optional
: 内部領域の境界. 値に応じてその格子点で境界値計算 1 = 固定端境界, 10 = 境界の内側. 2 = y 方向自由端境界 (フラックスは上向き) -2 = y 方向自由端境界 (フラックスは下向き) 4 = x 方向自由端境界 (フラックスは右向き) -4 = x 方向自由端境界 (フラックスは左向き) 3 = 周期境界 8 = |_, ~| で両方とも自由境界条件 -8 = |~, _| で両方とも自由境界条件 この引数が与えられなければ全領域を計算する. 境界の内側格子点 (10) は反復計算を行わず, undef で設定された値もしくはゼロが入る. このときの境界値は bound_opt の値が用いられる.
init_flag :logical, intent(in), optional
: psi の値をゼロで初期化するか. .true. = 初期化する. .false. = 初期化しない. デフォルトでは初期化する.
accel :real, intent(in), optional
: SOR の加速係数 (0 < accel < 2) デフォルト = 1
ln :integer, intent(in), optional
: 反復回数 この値が与えられるとき, eps の値に無関係に ln 回ループさせる.

ガウス=ザイデル法による楕円型方程式の求積 各オプション配列は, ポアソン系の各微分項の値. デフォルトはゼロで設定される. $$xa\dfrac{partial ^2\psi}{partial x^2} +ya\dfrac{partial ^2\psi}{partial y^2} +za\dfrac{partial ^2\psi}{partial z^2} +a\dfrac{partial ^2\psi}{partial x\partial y} +b\dfrac{partial ^2\psi}{partial y\partial z} +c\dfrac{partial ^2\psi}{partial z\partial x} +d\dfrac{partial psi}{partial x} +e\dfrac{partial psi}{partial y} +f\dfrac{partial psi}{partial z} +g\psi =rho $$ の各係数に対応している.

[Source]

subroutine Ellip_GauSei_3d(x, y, z, rho, eps, boundary, psi, bound_opt, xa, ya, za, a, b, c, d, e, f, g, undef, inner_bound, init_flag, accel, ln )
! ガウス=ザイデル法による楕円型方程式の求積
! 各オプション配列は, ポアソン系の各微分項の値. デフォルトはゼロで設定される.
! $$xa\dfrac{\partial ^2\psi}{\partial x^2} +ya\dfrac{\partial ^2\psi}{\partial y^2} +za\dfrac{\partial ^2\psi}{\partial z^2} +a\dfrac{\partial ^2\psi}{\partial x\partial y} +b\dfrac{\partial ^2\psi}{\partial y\partial z} +c\dfrac{\partial ^2\psi}{\partial z\partial x} +d\dfrac{\partial \psi}{\partial x} +e\dfrac{\partial \psi}{\partial y} +f\dfrac{\partial \psi}{\partial z} +g\psi =\rho $$
! の各係数に対応している.
  implicit none
  real, intent(in) :: x(:)  ! 領域の x 座標
  real, intent(in) :: y(:)  ! 領域の y 座標
  real, intent(in) :: z(:)  ! 領域の z 座標
  real, intent(in) :: rho(size(x),size(y),size(z))  ! ポアソン方程式の強制項
                   ! rho =0 でラプラス方程式も求積可能
  real, intent(in) :: eps  ! 収束条件
  character(6), intent(in) :: boundary  ! 境界条件
                ! 4 文字で各辺の境界条件を与える.
                ! 1 文字目 : x 下端, 2 文字目 : y 左端, 3 文字目 : x 上端,
                ! 4 文字目 : y 右端
                ! boundary は 1 : 固定端境界, 2 : 自由端境界, 3 : 周期境界
  real, intent(in), optional :: bound_opt(size(x),size(y),size(z))
                             ! 境界での強制
                             ! ノイマン境界の場合 : フラックス値
  real, intent(in), optional :: xa(size(x),size(y),size(z))  ! 各微分項の係数
  real, intent(in), optional :: ya(size(x),size(y),size(z))  ! 各微分項の係数
  real, intent(in), optional :: za(size(x),size(y),size(z))  ! 各微分項の係数
  real, intent(in), optional :: a(size(x),size(y),size(z))  ! 各微分項の係数
  real, intent(in), optional :: b(size(x),size(y),size(z))  ! 各微分項の係数
  real, intent(in), optional :: c(size(x),size(y),size(z))  ! 各微分項の係数
  real, intent(in), optional :: d(size(x),size(y),size(z))  ! 各微分項の係数
  real, intent(in), optional :: e(size(x),size(y),size(z))  ! 各微分項の係数
  real, intent(in), optional :: f(size(x),size(y),size(z))  ! 各微分項の係数
  real, intent(in), optional :: g(size(x),size(y),size(z))  ! 各微分項の係数
  real, intent(inout) :: psi(size(x),size(y),size(z))  ! ポアソン方程式の解
  real, intent(in), optional :: undef  ! 未定義値
  integer, intent(in), optional :: inner_bound(size(x),size(y),size(z))
                             ! 内部領域の境界. 値に応じてその格子点で境界値計算
                             ! 1 = 固定端境界, 10 = 境界の内側.
                             ! 2 = y 方向自由端境界 (フラックスは上向き)
                             ! -2 = y 方向自由端境界 (フラックスは下向き)
                             ! 4 = x 方向自由端境界 (フラックスは右向き)
                             ! -4 = x 方向自由端境界 (フラックスは左向き)
                             ! 3 = 周期境界
                             ! 8 = |_, ~| で両方とも自由境界条件
                             ! -8 = |~, _| で両方とも自由境界条件
                             ! この引数が与えられなければ全領域を計算する.
                             ! 境界の内側格子点 (10) は反復計算を行わず,
                             ! undef で設定された値もしくはゼロが入る.
                             ! このときの境界値は bound_opt の値が用いられる.
  logical, intent(in), optional :: init_flag  ! psi の値をゼロで初期化するか.
                             ! .true. = 初期化する. .false. = 初期化しない.
                             ! デフォルトでは初期化する.
  real, intent(in), optional :: accel  ! SOR の加速係数 (0 < accel < 2)
                             ! デフォルト = 1
  integer, intent(in), optional :: ln  ! 反復回数
                             ! この値が与えられるとき, eps の値に無関係に
                             ! ln 回ループさせる.

  integer :: i, j, k, ix, jy, kz, nl, counter
  integer :: nx  ! x 方向の配列要素
  integer :: ny  ! y 方向の配列要素
  integer :: nz  ! z 方向の配列要素
  integer :: signa, signb, signc  ! 各係数を計算するかどうか
  integer, dimension(size(x),size(y),size(z)) :: ib

  real :: defun
  real :: tmp, err, err_max
  real :: tmp_b, accc
  real :: bnd(size(x),size(y),size(z))
  real :: dx(size(x)), dy(size(y)), dx2(size(x)), dy2(size(y))
  real :: dz(size(z)), dz2(size(z))
  real, dimension(size(x),size(y)) :: dxdy
  real, dimension(size(y),size(z)) :: dydz
  real, dimension(size(x),size(z)) :: dxdz
  real, dimension(size(x),size(y),size(z)) :: xt, yt, zt, at, bt, ct, dt, et, ft, gt
  real, dimension(size(x),size(y),size(z)) :: xdp, xdm, yep, yem, zfp, zfm, xyz, divi

  character(6) :: bound
  logical :: sor_flag
  logical, dimension(size(x),size(y),size(z)) :: inner_flag

  bound(1:6)=boundary(1:6)

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

!-- 応答関数の初期化

  if(present(init_flag))then
     if(init_flag.eqv..true.)then
        psi = 0.0
     end if
  else
     psi = 0.0
  end if

!-- 内部境界の判別フラグの設定

  if(present(inner_bound))then
     call set_bound_3d( bound, ib, inner_flag, inner_bound )
  else
     call set_bound_3d( bound, ib, inner_flag )
  end if

!-- 領域・内部境界における境界値の設定

  if(present(bound_opt))then
     call setval_bound_3d( ib, bnd, psi, bound_opt )
  else
     call setval_bound_3d( ib, bnd, psi )
  end if

!-- 未定義値の設定

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

!-- 係数の代入
!-- xa, ya, za については, 値が入れられていなければ, 全配列に 1 を代入する.
  if(present(xa))then
     call set_coe_3d( xt, ext=xa )
  else
     call set_coe_3d( xt, def=1.0 )
  end if

  if(present(ya))then
     call set_coe_3d( yt, ext=ya )
  else
     call set_coe_3d( yt, def=1.0 )
  end if

  if(present(za))then
     call set_coe_3d( zt, ext=za )
  else
     call set_coe_3d( zt, def=1.0 )
  end if

  if(present(a))then
     call set_coe_3d( at, ext=a )
     signa=1
  else
     call set_coe_3d( at, def=0.0 )
     signa=0
  end if

  if(present(b))then
     call set_coe_3d( bt, ext=b )
     signb=1
  else
     call set_coe_3d( bt, def=0.0 )
     signb=0
  end if

  if(present(c))then
     signc=1
     call set_coe_3d( ct, ext=c )
  else
     call set_coe_3d( ct, def=0.0 )
     signc=0
  end if

  if(present(d))then
     call set_coe_3d( dt, ext=d )
  else
     call set_coe_3d( dt, def=0.0 )
  end if

  if(present(e))then
     call set_coe_3d( et, ext=e )
  else
     call set_coe_3d( et, def=0.0 )
  end if

  if(present(f))then
     call set_coe_3d( ft, ext=f )
  else
     call set_coe_3d( ft, def=0.0 )
  end if

  if(present(g))then
     call set_coe_3d( gt, ext=g )
  else
     call set_coe_3d( gt, def=0.0 )
  end if

!-- 最高階数における係数チェック. (係数がゼロでないか調べる.)

  call check_coe_3d( xt, 0.0 )
  call check_coe_3d( yt, 0.0 )
  call check_coe_3d( zt, 0.0 )

!-- 加速係数の判定

  if(present(accel))then
     accc=accel
     sor_flag=.true.
  else
     sor_flag=.false.
  end if

!-- ln の処理
  if(present(ln))then
     nl=ln
  else
     nl=0
  end if

!-- 以下で先に格子間隔等の 1 回計算でよいものを求めておく.
!-- これらは 1 方向のみで変化すればよい.
!-- 格子点間隔の計算
  do i=2,nx-1
     dx(i)=(x(i+1)-x(i-1))*0.5
     dx2(i)=dx(i)**2
  end do
  do j=2,ny-1
     dy(j)=(y(j+1)-y(j-1))*0.5
     dy2(j)=dy(j)**2
  end do
  do k=2,nz-1
     dz(k)=(z(k+1)-z(k-1))*0.5
     dz2(k)=dz(k)**2
  end do

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

  do j=1,ny
     do i=1,nx
        dxdy(i,j)=dx(i)*dy(j)
     end do
  end do
  do k=1,nz
     do j=1,ny
        dydz(j,k)=dy(j)*dz(k)
     end do
  end do
  do k=1,nz
     do i=1,nx
        dxdz(i,k)=dx(i)*dz(k)
     end do
  end do

!-- ポアソン係数の計算
!-- 実際にポアソン係数が使用されるのは, 境界を除いた 1 回り内側なので,
!-- 計算量削減のため, ループをそのようにしておく.

  xyz=0.0
  xdp=0.0
  xdm=0.0
  yep=0.0
  yem=0.0
  zfp=0.0
  zfm=0.0

!-- 最高次数係数 ac の計算
!$omp parallel default(shared)
!$omp do schedule(dynamic) private(i,j,k)

  do k=2,nz-1
     do j=2,ny-1
        do i=2,nx-1
           xyz(i,j,k)=1.0/(2.0*(xt(i,j,k)/dx2(i)+yt(i,j,k)/dy2(j)+zt(i,j,k)/dz2(k))-gt(i,j,k))
           xdp(i,j,k)=(xt(i,j,k)/(dx2(i))+0.5*dt(i,j,k)/dx(i))*xyz(i,j,k)
           xdm(i,j,k)=(xt(i,j,k)/(dx2(i))-0.5*dt(i,j,k)/dx(i))*xyz(i,j,k)
           yep(i,j,k)=(yt(i,j,k)/(dy2(j))+0.5*et(i,j,k)/dy(j))*xyz(i,j,k)
           yem(i,j,k)=(yt(i,j,k)/(dy2(j))-0.5*et(i,j,k)/dy(j))*xyz(i,j,k)
           zfp(i,j,k)=(zt(i,j,k)/(dz2(k))+0.5*ft(i,j,k)/dz(k))*xyz(i,j,k)
           zfm(i,j,k)=(zt(i,j,k)/(dz2(k))-0.5*ft(i,j,k)/dz(k))*xyz(i,j,k)
           if(signa==1)then
              at(i,j,k)=0.25*at(i,j,k)/(dxdy(i,j))*xyz(i,j,k)
           else
              at(i,j,k)=0.0
           end if
           if(signb==1)then
              bt(i,j,k)=0.25*bt(i,j,k)/(dydz(j,k))*xyz(i,j,k)
           else
              bt(i,j,k)=0.0
           end if
           if(signc==1)then
              ct(i,j,k)=0.25*ct(i,j,k)/(dxdz(i,k))*xyz(i,j,k)
           else
              ct(i,j,k)=0.0
           end if
        end do
     end do
  end do

!$omp end do
!$omp end parallel

  err_max=eps  ! while に入るための便宜的措置
  counter=0

!-- 実際のソルバ ---
  do while(err_max>=eps)
     err_max=0.0

     do k=1,nz
        do j=1,ny
           do i=1,nx

!-- 以降, 反復計算に必要な点の周囲 8 点についてそれぞれ
!-- inner_flag のチェックと同時に適切な値をそれぞれ逐次計算する,
!-- 下で, 各方向の格子計算を行っているが, ib=8,-8 は 4 隅格子のときにしか
!-- case select しない. なぜなら, 真横や上下に隅格子がくることはありえない.
              if(inner_flag(i,j,k).eqv..false.)then  ! .false. なら領域計算開始

                 tmp=-rho(i,j,k)*xyz(i,j,k)
                 tmp=tmp+xdp(i,j,k)*psi(i+1,j,k) +xdm(i,j,k)*psi(i-1,j,k) +yep(i,j,k)*psi(i,j+1,k) +yem(i,j,k)*psi(i,j-1,k) +zfp(i,j,k)*psi(i,j,k+1) +zfm(i,j,k)*psi(i,j,k-1)

                 if(signa==1)then
                    tmp_b=at(i,j,k)*(psi(i+1,j+1,k)+psi(i-1,j-1,k) -psi(i+1,j-1,k)-psi(i-1,j+1,k))
                 else
                    tmp_b=0.0
                 end if
                 if(signb==1)then  ! そもそも bt = 0 なら計算しない.
                    tmp_b=tmp_b+bt(i,j,k)*(psi(i,j+1,k+1)+psi(i,j-1,k-1) -psi(i,j-1,k+1)-psi(i,j+1,k-1))
                 end if
                 if(signc==1)then  ! そもそも bt = 0 なら計算しない.
                    tmp_b=tmp_b+ct(i,j,k)*(psi(i+1,j,k+1)+psi(i-1,j,k-1) -psi(i-1,j,k+1)-psi(i+1,j,k-1))
                 end if
              
                 tmp=tmp+tmp_b

              else   ! .true. なら境界計算開始.

                 select case (ib(i,j,k))
                 case (1)
                    tmp=bnd(i,j,k)

                 case (2)  ! x 方向にフラックス一定, 上側が参照値
                    tmp=psi(i+1,j,k)-bnd(i,j,k)*dx(i)

                 case (-2)  ! x 方向にフラックス一定, 下側が参照値
                    tmp=psi(i-1,j,k)+bnd(i,j,k)*dx(i)

                 case (4)  ! y 方向にフラックス一定, 右側が参照値
                    tmp=psi(i,j+1,k)-bnd(i,j,k)*dy(j)

                 case (-4)  ! y 方向にフラックス一定, 左側が参照値
                    tmp=psi(i,j-1,k)+bnd(i,j,k)*dy(j)

                 case (6)  ! y 方向にフラックス一定, 右側が参照値
                    tmp=psi(i,j,k+1)-bnd(i,j,k)*dz(k)

                 case (-6)  ! y 方向にフラックス一定, 左側が参照値
                    tmp=psi(i,j,k-1)+bnd(i,j,k)*dz(k)

                 case (3)  ! 12 辺, もしくは 8 点で周期境界を判断
                    if(i==1)then
                       ix=nx-1
                    else if(i==nx)then
                       ix=2
                    else
                       ix=i
                    end if
                    if(j==1)then
                       jy=ny-1
                    else if(j==ny)then
                       jy=2
                    else
                       jy=j
                    end if
                    if(k==1)then
                       kz=nz-1
                    else if(k==nz)then
                       kz=2
                    else
                       kz=k
                    end if
                    tmp=psi(ix,jy,kz)

                 case (8)  ! 両方フラックス一定で z 面の x, y 右上か左下角.
                    if(i==1.and.j==1)then  ! -- 評価 1
                       tmp=psi(i+1,j+1,k) -0.5*(bnd(i+1,j,k)*dy(j)+bnd(i,j+1,k)*dx(i))

                    else if(i==nx.and.j==ny)then  ! -- 評価 2
                       tmp=psi(i-1,j-1,k) +0.5*(bnd(i-1,j,k)*dy(j)+bnd(i,j-1,k)*dx(i))

                    else if(ib(i-1,j,k)==10.and.ib(i,j-1,k)==10)then
                       ! -- 評価 1 と同じ
                       tmp=psi(i+1,j+1,k) -0.5*(bnd(i+1,j,k)*dy(j)+bnd(i,j+1,k)*dx(i))

                    else if(ib(i+1,j,k)==10.and.ib(i,j+1,k)==10)then
                       ! -- 評価 2 と同じ
                       tmp=psi(i-1,j-1,k) +0.5*(bnd(i-1,j,k)*dy(j)+bnd(i,j-1,k)*dx(i))

                    end if

                 case (-8)  ! 両方フラックス一定で z 面の x, y 右下か左上角.
                    if(i==1.and.j==ny)then  ! -- 評価 1
                       tmp=psi(i+1,j-1,k) +0.5*(bnd(i+1,j,k)*dy(j)-bnd(i,j-1,k)*dx(i))

                    else if(i==nx.and.j==1)then  ! -- 評価 2
                       tmp=psi(i-1,j+1,k) +0.5*(-bnd(i-1,j,k)*dy(j)+bnd(i,j+1,k)*dx(i))

                    else if(ib(i-1,j,k)==10.and.ib(i,j+1,k)==10)then
                       ! -- 評価 1 と同じ
                       tmp=psi(i+1,j-1,k) +0.5*(bnd(i+1,j,k)*dy(j)-bnd(i,j-1,k)*dx(i))

                    else if(ib(i+1,j,k)==10.and.ib(i,j-1,k)==10)then
                       ! -- 評価 2 と同じ
                       tmp=psi(i-1,j+1,k) +0.5*(-bnd(i-1,j,k)*dy(j)+bnd(i,j+1,k)*dx(i))
                    end if

                 case (12)  ! 両方フラックス一定で y 面の x, z 右上か左下角.
                    if(i==1.and.k==1)then  ! -- 評価 1
                       tmp=psi(i+1,j,k+1) -0.5*(bnd(i+1,j,k)*dz(k)+bnd(i,j,k+1)*dx(i))

                    else if(i==nx.and.k==nz)then  ! -- 評価 2
                       tmp=psi(i-1,j,k-1) +0.5*(bnd(i-1,j,k)*dz(k)+bnd(i,j,k-1)*dx(i))

                    else if(ib(i-1,j,k)==10.and.ib(i,j,k-1)==10)then
                       ! -- 評価 1 と同じ
                       tmp=psi(i+1,j,k+1) -0.5*(bnd(i+1,j,k)*dz(k)+bnd(i,j,k+1)*dx(i))

                    else if(ib(i+1,j,k)==10.and.ib(i,j,k+1)==10)then
                       ! -- 評価 2 と同じ
                       tmp=psi(i-1,j,k-1) +0.5*(bnd(i-1,j,k)*dz(k)+bnd(i,j,k-1)*dx(i))
                    end if

                 case (-12)  ! 両方フラックス一定で y 面の x, z 右下か左上角.
                    if(i==1.and.k==nz)then  ! -- 評価 1
                       tmp=psi(i+1,j,k-1) +0.5*(bnd(i+1,j,k)*dz(k)-bnd(i,j,k-1)*dx(i))

                    else if(i==nx.and.k==1)then  ! -- 評価 2
                       tmp=psi(i-1,j,k+1) +0.5*(-bnd(i-1,j,k)*dz(k)+bnd(i,j,k+1)*dx(i))

                    else if(ib(i-1,j,k)==10.and.ib(i,j,k+1)==10)then
                       ! -- 評価 1 と同じ
                       tmp=psi(i+1,j,k-1) +0.5*(bnd(i+1,j,k)*dz(k)-bnd(i,j,k-1)*dx(i))

                    else if(ib(i+1,j,k)==10.and.ib(i,j,k-1)==10)then
                       ! -- 評価 2 と同じ
                       tmp=psi(i-1,j,k+1) +0.5*(-bnd(i-1,j,k)*dz(k)+bnd(i,j,k+1)*dx(i))
                    end if

                 case (24)  ! 両方フラックス一定で x 面の y, z 右上か左下角.
                    if(j==1.and.k==1)then  ! -- 評価 1
                       tmp=psi(i,j+1,k+1) -0.5*(bnd(i,j+1,k)*dz(k)+bnd(i,j,k+1)*dy(j))

                    else if(j==ny.and.k==nz)then  ! -- 評価 2
                       tmp=psi(i,j-1,k-1) +0.5*(bnd(i,j-1,k)*dz(k)+bnd(i,j,k-1)*dy(j))

                    else if(ib(i,j-1,k)==10.and.ib(i,j,k-1)==10)then
                       ! -- 評価 1 と同じ
                       tmp=psi(i,j+1,k+1) -0.5*(bnd(i,j+1,k)*dz(k)+bnd(i,j,k+1)*dy(j))

                    else if(ib(i,j+1,k)==10.and.ib(i,j,k+1)==10)then
                       ! -- 評価 2 と同じ
                       tmp=psi(i,j-1,k-1) +0.5*(bnd(i,j-1,k)*dz(k)+bnd(i,j,k-1)*dy(j))
                    end if

                 case (-24)  ! 両方フラックス一定で x 面の y, z 右下か左上角.
                    if(j==1.and.k==nz)then  ! -- 評価 1
                       tmp=psi(i,j+1,k-1) +0.5*(bnd(i,j+1,k)*dz(k)-bnd(i,j,k-1)*dy(j))

                    else if(j==ny.and.k==1)then  ! -- 評価 2
                       tmp=psi(i,j-1,k+1) +0.5*(-bnd(i,j-1,k)*dz(k)+bnd(i,j,k+1)*dy(j))

                    else if(ib(i,j-1,k)==10.and.ib(i,j,k+1)==10)then
                       ! -- 評価 1 と同じ
                       tmp=psi(i,j+1,k-1) +0.5*(bnd(i,j+1,k)*dz(k)-bnd(i,j,k-1)*dy(j))

                    else if(ib(i,j+1,k)==10.and.ib(i,j,k-1)==10)then
                       ! -- 評価 2 と同じ
                       tmp=psi(i,j-1,k+1) +0.5*(-bnd(i,j-1,k)*dz(k)+bnd(i,j,k+1)*dy(j))
                    end if

                 !-- 以降, 隅領域なので, 個別に設定.
                 case (11)  ! 両方フラックス一定で (1,1,1) 点, もしくは内部領域.
                    tmp=psi(i+1,j+1,k+1) -(bnd(i,j+1,k+1)*dx(i) +bnd(i+1,j,k+1)*dy(j) +bnd(i+1,j+1,k)*dz(k))/3.0

                 case (13)  ! 両方フラックス一定で (nx,1,1) 点, もしくは内部領域.
                    tmp=psi(i-1,j+1,k+1) -(-bnd(i,j+1,k+1)*dx(i) +bnd(i-1,j,k+1)*dy(j) +bnd(i-1,j+1,k)*dz(k))/3.0

                 case (17)  ! 両方フラックス一定で (1,ny,1) 点, もしくは内部領域.
                    tmp=psi(i+1,j-1,k+1) -(bnd(i,j-1,k+1)*dx(i) -bnd(i+1,j,k+1)*dy(j) +bnd(i+1,j-1,k)*dz(k))/3.0

                 case (19)  ! 両方フラックス一定で (nx,ny,1) 点, もしくは内部領域.
                    tmp=psi(i-1,j-1,k+1) -(-bnd(i,j-1,k+1)*dx(i) -bnd(i-1,j,k+1)*dy(j) +bnd(i-1,j-1,k)*dz(k))/3.0

                 case (23)  ! 両方フラックス一定で (1,1,nz) 点, もしくは内部領域.
                    tmp=psi(i+1,j+1,k-1) -(bnd(i,j+1,k-1)*dx(i) +bnd(i+1,j,k-1)*dy(j) -bnd(i+1,j+1,k)*dz(k))/3.0

                 case (29)  ! 両方フラックス一定で (nx,1,nz) 点, もしくは内部領域.
                    tmp=psi(i-1,j+1,k-1) -(-bnd(i,j+1,k-1)*dx(i) +bnd(i-1,j,k-1)*dy(j) -bnd(i-1,j+1,k)*dz(k))/3.0

                 case (31)  ! 両方フラックス一定で (1,ny,nz) 点, もしくは内部領域.
                    tmp=psi(i+1,j-1,k-1) -(bnd(i,j-1,k-1)*dx(i) -bnd(i+1,j,k-1)*dy(j) -bnd(i+1,j-1,k)*dz(k))/3.0

                 case (37)  ! 両方フラックス一定で (nx,ny,nz) 点, もしくは内部領域.
                    tmp=psi(i-1,j-1,k-1) -(-bnd(i,j-1,k-1)*dx(i) -bnd(i-1,j,k-1)*dy(j) -bnd(i-1,j-1,k)*dz(k))/3.0

                 case (-11)  ! 両方フラックス一定で (i+1,j+1,k+1) に 10 が設定.
                    tmp=(psi(i-1,j+1,k+1)+psi(i+1,j-1,k+1)+psi(i+1,j+1,k-1) +bnd(i,j+1,k+1)*dx(i)+bnd(i+1,j,k+1)*dy(j) +bnd(i+1,j+1,k)*dz(k))/3.0

                 case (-13)  ! 両方フラックス一定で (i-1,j+1,k+1) に 10 が設定.
                    tmp=(psi(i+1,j+1,k+1)+psi(i-1,j-1,k+1)+psi(i-1,j+1,k-1) -bnd(i,j+1,k+1)*dx(i)+bnd(i-1,j,k+1)*dy(j) +bnd(i-1,j+1,k)*dz(k))/3.0

                 case (-17)  ! 両方フラックス一定で (i+1,j-1,k+1) に 10 が設定.
                    tmp=(psi(i-1,j-1,k+1)+psi(i+1,j+1,k+1)+psi(i+1,j-1,k-1) +bnd(i,j-1,k+1)*dx(i)-bnd(i+1,j,k+1)*dy(j) +bnd(i+1,j-1,k)*dz(k))/3.0

                 case (-19)  ! 両方フラックス一定で (i-1,j-1,k+1) に 10 が設定.
                    tmp=(psi(i+1,j-1,k+1)+psi(i-1,j+1,k+1)+psi(i-1,j-1,k-1) -bnd(i,j-1,k+1)*dx(i)-bnd(i-1,j,k+1)*dy(j) +bnd(i-1,j-1,k)*dz(k))/3.0

                 case (-23)  ! 両方フラックス一定で (i+1,j+1,k-1) に 10 が設定.
                    tmp=(psi(i-1,j+1,k-1)+psi(i+1,j-1,k-1)+psi(i+1,j+1,k+1) +bnd(i,j+1,k-1)*dx(i)+bnd(i+1,j,k-1)*dy(j) -bnd(i+1,j+1,k)*dz(k))/3.0

                 case (-29)  ! 両方フラックス一定で (i-1,j+1,k-1) に 10 が設定.
                    tmp=(psi(i+1,j+1,k-1)+psi(i-1,j-1,k-1)+psi(i-1,j+1,k+1) -bnd(i,j+1,k-1)*dx(i)+bnd(i-1,j,k-1)*dy(j) -bnd(i-1,j+1,k)*dz(k))/3.0

                 case (-31)  ! 両方フラックス一定で (i+1,j-1,k-1) に 10 が設定.
                    tmp=(psi(i-1,j-1,k-1)+psi(i+1,j+1,k-1)+psi(i+1,j-1,k+1) +bnd(i,j-1,k-1)*dx(i)-bnd(i+1,j,k-1)*dy(j) -bnd(i+1,j-1,k)*dz(k))/3.0

                 case (-37)  ! 両方フラックス一定で (i-1,j-1,k-1) に 10 が設定.
                    tmp=(psi(i+1,j-1,k-1)+psi(i-1,j+1,k-1)+psi(i-1,j-1,k+1) -bnd(i,j-1,k-1)*dx(i)-bnd(i-1,j,k-1)*dy(j) -bnd(i-1,j-1,k)*dz(k))/3.0

                 end select

              end if

              if(sor_flag.eqv..true.)then
                 tmp=(1.0-accc)*psi(i,j,k)+tmp*accc
              end if

              err=abs(tmp-psi(i,j,k))

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

              psi(i,j,k)=tmp

           end do
        end do
     end do

     if(nl/=0)then
        counter=counter+1
        if(counter==nl)then
           exit
        else
           err_max=eps
        end if
     end if

  end do

!-- 境界の設定

  call calculate_bound_3d( ib, dx, dy, dz, bnd, psi )

!-- 未定義領域には undef を代入する.

  do k=1,nz
     do j=1,ny
        do i=1,nx
           if(ib(i,j,k)==10)then
              psi(i,j,k)=defun
           end if
        end do
     end do
  end do

end subroutine Ellip_GauSei_3d
Subroutine :
x(:) :real, intent(in)
: 領域の横座標
y(:) :real, intent(in)
: 領域の縦座標
rho(size(x),size(y)) :real, intent(in)
: ポアソン方程式の強制項 rho =0 でラプラス方程式も求積可能
eps :real, intent(in)
: 収束条件
boundary :character(4), intent(in)
: 境界条件 4 文字で各辺の境界条件を与える. 1 文字目 : x 下端, 2 文字目 : y 左端, 3 文字目 : x 上端, 4 文字目 : y 右端 boundary は 1 : 固定端境界, 2 : 自由端境界, 3 : 周期境界
psi(size(x),size(y)) :real, intent(inout)
: ポアソン方程式の解
bound_opt(size(x),size(y)) :real, intent(in), optional
: 境界での強制 ノイマン境界の場合 : フラックス値
a(size(x),size(y)) :real, intent(in), optional
: 各微分項の係数
b(size(x),size(y)) :real, intent(in), optional
: 各微分項の係数
c(size(x),size(y)) :real, intent(in), optional
: 各微分項の係数
d(size(x),size(y)) :real, intent(in), optional
: 各微分項の係数
e(size(x),size(y)) :real, intent(in), optional
: 各微分項の係数
f(size(x),size(y)) :real, intent(in), optional
: 各微分項の係数
undef :real, intent(in), optional
: 未定義値
inner_bound(size(x),size(y)) :integer, intent(in), optional
: 内部領域の境界. 値に応じてその格子点で境界値計算 1 = 固定端境界, 10 = 境界の内側. 2 = y 方向自由端境界 (フラックスは上向き) -2 = y 方向自由端境界 (フラックスは下向き) 4 = x 方向自由端境界 (フラックスは右向き) -4 = x 方向自由端境界 (フラックスは左向き) 3 = 周期境界, -3 = 隅領域で両方とも周期境界 8 = |_, ~| で両方とも自由境界条件 -8 = |~, _| で両方とも自由境界条件 この引数が与えられなければ全領域を計算する. 境界の内側格子点 (10) は反復計算を行わず, undef で設定された値もしくはゼロが入る. このときの境界値は bound_opt の値が用いられる.
init_flag :logical, intent(in), optional
: psi の値をゼロで初期化するか. .true. = 初期化する. .false. = 初期化しない. デフォルトでは初期化する.
accel :real, intent(in), optional
: SOR の加速係数 (0 < accel < 2) デフォルト = 1
ln :integer, intent(in), optional
: 反復回数 この値が与えられるとき, eps の値に無関係に ln 回ループさせる.

openmp によるスレッド並列が可能. ガウス・ザイデル法ではアルゴリズムの点から並列化が困難と思われたので, 並列計算によるポアソン方程式の求積が必要となるなら, ヤコビ法のものを使用されたい.

[Source]

subroutine Ellip_Jacobi_2d( x, y, rho, eps, boundary, psi, bound_opt, a, b, c, d, e, f, undef, inner_bound, init_flag, accel, ln )
  ! openmp によるスレッド並列が可能.
  ! ガウス・ザイデル法ではアルゴリズムの点から並列化が困難と思われたので,
  ! 並列計算によるポアソン方程式の求積が必要となるなら,
  ! ヤコビ法のものを使用されたい.
  implicit none
  real, intent(in) :: x(:)  ! 領域の横座標
  real, intent(in) :: y(:)  ! 領域の縦座標
  real, intent(in) :: rho(size(x),size(y))  ! ポアソン方程式の強制項
                   ! rho =0 でラプラス方程式も求積可能
  real, intent(in) :: eps  ! 収束条件
  character(4), intent(in) :: boundary  ! 境界条件
                ! 4 文字で各辺の境界条件を与える.
                ! 1 文字目 : x 下端, 2 文字目 : y 左端, 3 文字目 : x 上端,
                ! 4 文字目 : y 右端
                ! boundary は 1 : 固定端境界, 2 : 自由端境界, 3 : 周期境界
  real, intent(in), optional :: bound_opt(size(x),size(y))  ! 境界での強制
                             ! ノイマン境界の場合 : フラックス値
  real, intent(in), optional :: a(size(x),size(y))  ! 各微分項の係数
  real, intent(in), optional :: b(size(x),size(y))  ! 各微分項の係数
  real, intent(in), optional :: c(size(x),size(y))  ! 各微分項の係数
  real, intent(in), optional :: d(size(x),size(y))  ! 各微分項の係数
  real, intent(in), optional :: e(size(x),size(y))  ! 各微分項の係数
  real, intent(in), optional :: f(size(x),size(y))  ! 各微分項の係数
  real, intent(in), optional :: undef  ! 未定義値
  integer, intent(in), optional :: inner_bound(size(x),size(y))
                             ! 内部領域の境界. 値に応じてその格子点で境界値計算
                             ! 1 = 固定端境界, 10 = 境界の内側.
                             ! 2 = y 方向自由端境界 (フラックスは上向き)
                             ! -2 = y 方向自由端境界 (フラックスは下向き)
                             ! 4 = x 方向自由端境界 (フラックスは右向き)
                             ! -4 = x 方向自由端境界 (フラックスは左向き)
                             ! 3 = 周期境界, -3 = 隅領域で両方とも周期境界
                             ! 8 = |_, ~| で両方とも自由境界条件
                             ! -8 = |~, _| で両方とも自由境界条件
                             ! この引数が与えられなければ全領域を計算する.
                             ! 境界の内側格子点 (10) は反復計算を行わず,
                             ! undef で設定された値もしくはゼロが入る.
                             ! このときの境界値は bound_opt の値が用いられる.
  real, intent(inout) :: psi(size(x),size(y))  ! ポアソン方程式の解
  logical, intent(in), optional :: init_flag  ! psi の値をゼロで初期化するか.
                             ! .true. = 初期化する. .false. = 初期化しない.
                             ! デフォルトでは初期化する.
  real, intent(in), optional :: accel  ! SOR の加速係数 (0 < accel < 2)
                             ! デフォルト = 1
  integer, intent(in), optional :: ln  ! 反復回数
                             ! この値が与えられるとき, eps の値に無関係に
                             ! ln 回ループさせる.

  integer :: i, j, ix, jy, nl, counter
  integer :: nx  ! x 方向の配列要素
  integer :: ny  ! y 方向の配列要素
  integer :: signb  ! クロスターム b が存在するか
  integer, dimension(size(x),size(y)) :: ib

  real :: defun
  real :: err, err_max, accc
  real :: bnd(size(x),size(y))
  real :: dx(size(x)), dy(size(y)), dx2(size(x)), dy2(size(y))
  real, dimension(size(x),size(y)) :: dxdy
  real, dimension(size(x),size(y)) :: at, bt, ct, dt, et, ft
  real, dimension(size(x),size(y)) :: adp, adm, cep, cem, ac
  real, dimension(size(x),size(y)) :: tmp, tmp_b, divi

  character(4) :: bound
  logical :: sor_flag
  logical, dimension(size(x),size(y)) :: inner_flag

  bound(1:4)=boundary(1:4)

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

!-- 応答関数の初期化

  if(present(init_flag))then
     if(init_flag.eqv..true.)then
        psi = 0.0
     end if
  else
     psi = 0.0
  end if

!-- 内部境界の判別フラグの設定

  if(present(inner_bound))then
     call set_bound( bound, ib, inner_flag, inner_bound )
  else
     call set_bound( bound, ib, inner_flag )
  end if

!-- 領域・内部境界における境界値の設定

  if(present(bound_opt))then
     call setval_bound( ib, bnd, psi, bound_opt )
  else
     call setval_bound( ib, bnd, psi )
  end if

!-- 未定義値の設定

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

!-- 係数の代入
!-- a, c については, 値が入れられていなければ, 全配列に 1 を代入する.

  if(present(a))then
     call set_coe( at, ext=a )
  else
     call set_coe( at, def=1.0 )
  end if

  if(present(c))then
     call set_coe( ct, ext=c )
  else
     call set_coe( ct, def=1.0 )
  end if

  if(present(b))then
     call set_coe( bt, ext=b )
     signb=1
  else
     call set_coe( bt, def=0.0 )
     signb=0
  end if

  if(present(d))then
     call set_coe( dt, ext=d )
  else
     call set_coe( dt, def=0.0 )
  end if

  if(present(e))then
     call set_coe( et, ext=e )
  else
     call set_coe( et, def=0.0 )
  end if

  if(present(f))then
     call set_coe( ft, ext=f )
  else
     call set_coe( ft, def=0.0 )
  end if

!-- 最高階数における係数チェック. (係数がゼロでないか調べる.)

  call check_coe( at, 0.0 )
  call check_coe( ct, 0.0 )

!-- 加速係数の判定

  if(present(accel))then
     accc=accel
     sor_flag=.true.
  else
     sor_flag=.false.
  end if

!-- ln の処理
  if(present(ln))then
     nl=ln
  else
     nl=0
  end if

!-- 以下で先に格子間隔等の 1 回計算でよいものを求めておく.
!-- これらは 1 方向のみで変化すればよい.
!-- 格子点間隔の計算
  do i=2,nx-1
     dx(i)=(x(i+1)-x(i-1))*0.5
     dx2(i)=dx(i)**2
  end do
  do j=2,ny-1
     dy(j)=(y(j+1)-y(j-1))*0.5
     dy2(j)=dy(j)**2
  end do

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

  do j=1,ny
     do i=1,nx
        dxdy(i,j)=dx(i)*dy(j)
     end do
  end do

!-- ポアソン係数の計算
!-- 実際にポアソン係数が使用されるのは, 境界を除いた 1 回り内側なので,
!-- 計算量削減のため, ループをそのようにしておく.

  ac=0.0
  adp=0.0
  adm=0.0
  cep=0.0
  cem=0.0
  bt=0.0

!-- 最高次数係数 ac の計算
!$omp parallel default(shared)
!$omp do schedule(dynamic) private(i,j)

  do j=2,ny-1
     do i=2,nx-1
        ac(i,j)=1.0/(2.0*(at(i,j)/dx2(i)+ct(i,j)/dy2(j))-ft(i,j))
        adp(i,j)=(at(i,j)/(dx2(i))+0.5*dt(i,j)/dx(i))*ac(i,j)
        adm(i,j)=(at(i,j)/(dx2(i))-0.5*dt(i,j)/dx(i))*ac(i,j)
        cep(i,j)=(ct(i,j)/(dy2(j))+0.5*et(i,j)/dy(j))*ac(i,j)
        cem(i,j)=(ct(i,j)/(dy2(j))-0.5*et(i,j)/dy(j))*ac(i,j)
        bt(i,j)=0.25*bt(i,j)/(dxdy(i,j))*ac(i,j)
     end do
  end do

!$omp end do
!$omp end parallel

  err_max=eps  ! while に入るための便宜的措置
  counter=0

!-- 実際のソルバ ---
  do while(err_max>=eps)
     err_max=0.0
!$omp parallel default(shared)
!$omp do schedule(dynamic) private(i,j,ix,jy)
     do j=1,ny
        do i=1,nx

!-- 以降, 反復計算に必要な点の周囲 8 点についてそれぞれ
!-- inner_flag のチェックと同時に適切な値をそれぞれ逐次計算する,
!-- 下で, 各方向の格子計算を行っているが, ib=8,-8 は 4 隅格子のときにしか
!-- case select しない. なぜなら, 真横や上下に隅格子がくることはありえない.
           if(inner_flag(i,j).eqv..false.)then  ! .false. なら領域計算開始
              tmp(i,j)=-rho(i,j)*ac(i,j)
              tmp(i,j)=tmp(i,j)+adp(i,j)*psi(i+1,j) +adm(i,j)*psi(i-1,j) +cep(i,j)*psi(i,j+1) +cem(i,j)*psi(i,j-1)
              if(signb==0)then  ! そもそも bt = 0 なら計算しない.
                  tmp_b(i,j)=0.0
              else
                  tmp_b(i,j)=bt(i,j)*(psi(i+1,j+1)+psi(i-1,j-1) -psi(i+1,j-1)-psi(i-1,j+1))
              end if

              tmp(i,j)=tmp(i,j)+tmp_b(i,j)

           else  ! .true. なら境界計算に移行.

              select case (ib(i,j))
              case (1)
                 tmp(i,j)=bnd(i,j)

              case (2)  ! x 方向にフラックス一定, 上側が参照値
                 tmp(i,j)=psi(i+1,j)-bnd(i,j)*dx(i)

              case (-2)  ! x 方向にフラックス一定, 下側が参照値
                 tmp(i,j)=psi(i-1,j)+bnd(i,j)*dx(i)

              case (4)  ! y 方向にフラックス一定, 右側が参照値
                 tmp(i,j)=psi(i,j+1)-bnd(i,j)*dy(j)

              case (-4)  ! y 方向にフラックス一定, 左側が参照値
                 tmp(i,j)=psi(i,j-1)+bnd(i,j)*dy(j)

              case (3)  ! 周期境界
                 if(i==1)then
                    ix=nx-1
                 else if(i==nx)then
                    ix=2
                 else
                    ix=i
                 end if
                 if(j==1)then
                    jy=ny-1
                 else if(j==ny)then
                    jy=2
                 else
                    jy=j
                 end if

                 tmp(i,j)=psi(ix,jy)

              case (7)  ! 両方フラックス一定で内部境界限定.
                 if((ib(i+1,j+1)==10).and.(ib(i+1,j)/=10).and. (ib(i,j+1)/=10))then
                    tmp(i,j)=0.5*(psi(i-1,j)+psi(i,j-1)) +0.5*bnd(i,j)*(dy(j)+dx(i))

                 else if((ib(i-1,j-1)==10).and.(ib(i-1,j)/=10).and. (ib(i,j-1)/=10))then
                    tmp(i,j)=0.5*(psi(i+1,j)+psi(i,j+1)) -0.5*bnd(i,j)*(dy(j)+dx(i))

                 end if

              case (-7)  ! 両方フラックス一定で内部境界限定.
                 if((ib(i-1,j+1)==10).and.(ib(i-1,j)/=10).and. (ib(i,j+1)/=10))then
                    tmp(i,j)=0.5*(psi(i+1,j)+psi(i,j-1)) +0.5*bnd(i,j)*(dy(j)-dx(i))

                 else if((ib(i+1,j-1)==10).and.(ib(i+1,j)/=10).and. (ib(i,j-1)/=10))then
                    tmp(i,j)=0.5*(psi(i-1,j)+psi(i,j+1)) +0.5*bnd(i,j)*(-dy(j)+dx(i))

                 end if

              case (8)  ! 両方フラックス一定で左下角か右上角, もしくは内部境界.
                 if(i==1.and.j==1)then  ! -- 評価 1
                    tmp(i,j)=psi(i+1,j+1)-0.5*(bnd(i+1,j)*dy(j)+bnd(i,j+1)*dx(i))

                 else if(i==nx.and.j==ny)then  ! -- 評価 2
                    tmp(i,j)=psi(i-1,j-1)+0.5*(bnd(i-1,j)*dy(j)+bnd(i,j-1)*dx(i))

                 else if(ib(i-1,j)==10.and.ib(i,j-1)==10)then
                    ! -- 評価 1 と同じ
                    tmp(i,j)=psi(i+1,j+1)-0.5*(bnd(i+1,j)*dy(j)+bnd(i,j+1)*dx(i))

                 else if(ib(i+1,j)==10.and.ib(i,j+1)==10)then
                    ! -- 評価 2 と同じ
                    tmp(i,j)=psi(i-1,j-1)+0.5*(bnd(i-1,j)*dy(j)+bnd(i,j-1)*dx(i))

                 end if

              case (-8)  ! 両方フラックス一定で右下角か左上角
                 if(i==1.and.j==ny)then  ! -- 評価 1
                    tmp(i,j)=psi(i+1,j-1)+0.5*(bnd(i+1,j)*dy(j)-bnd(i,j-1)*dx(i))

                 else if(i==nx.and.j==1)then  ! -- 評価 2
                    tmp(i,j)=psi(i-1,j+1)+0.5*(-bnd(i-1,j)*dy(j)+bnd(i,j+1)*dx(i))

                 else if(ib(i-1,j)==10.and.ib(i,j+1)==10)then
                    ! -- 評価 1 と同じ
                    tmp(i,j)=psi(i+1,j-1)+0.5*(bnd(i+1,j)*dy(j)-bnd(i,j-1)*dx(i))

                 else if(ib(i+1,j)==10.and.ib(i,j-1)==10)then
                    ! -- 評価 2 と同じ
                    tmp(i,j)=psi(i-1,j+1)+0.5*(-bnd(i-1,j)*dy(j)+bnd(i,j+1)*dx(i))

                 end if
              end select

           end if

           if(sor_flag.eqv..true.)then
              tmp(i,j)=(1.0-accc)*psi(i,j)+tmp(i,j)*accc
           end if

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

!!-- ここまでは if 文のインデントを調節しない.
!        end if
!
!        tmp=tmp/ac(i,j)
!-- 誤差の計算 ---
     do j=2,ny-1
        do i=2,nx-1
           err=abs(tmp(i,j)-psi(i,j))

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

!-- 一斉更新

     do j=1,ny
        do i=1,nx
           psi(i,j)=tmp(i,j)
        end do
     end do

     if(nl/=0)then
        counter=counter+1
        if(counter==nl)then
           exit
        else
           err_max=eps
        end if
     end if

  end do
  
!-- 境界の設定

  call calculate_bound( ib, dx, dy, bnd, psi )

!-- 未定義領域には undef を代入する.

  do j=1,ny
     do i=1,nx
        if(ib(i,j)==10)then
           psi(i,j)=defun
        end if
     end do
  end do

end subroutine Ellip_Jacobi_2d
Subroutine :
x(:) :real, intent(in)
: 領域の x 座標
y(:) :real, intent(in)
: 領域の y 座標
z(:) :real, intent(in)
: 領域の z 座標
rho(size(x),size(y),size(z)) :real, intent(in)
: ポアソン方程式の強制項 rho =0 でラプラス方程式も求積可能
eps :real, intent(in)
: 収束条件
boundary :character(6), intent(in)
: 境界条件 4 文字で各辺の境界条件を与える. 1 文字目 : x 下端, 2 文字目 : y 左端, 3 文字目 : x 上端, 4 文字目 : y 右端 boundary は 1 : 固定端境界, 2 : 自由端境界, 3 : 周期境界
psi(size(x),size(y),size(z)) :real, intent(inout)
: ポアソン方程式の解
bound_opt(size(x),size(y),size(z)) :real, intent(in), optional
: 境界での強制 ノイマン境界の場合 : フラックス値
xa(size(x),size(y),size(z)) :real, intent(in), optional
: 各微分項の係数
ya(size(x),size(y),size(z)) :real, intent(in), optional
: 各微分項の係数
za(size(x),size(y),size(z)) :real, intent(in), optional
: 各微分項の係数
a(size(x),size(y),size(z)) :real, intent(in), optional
: 各微分項の係数
b(size(x),size(y),size(z)) :real, intent(in), optional
: 各微分項の係数
c(size(x),size(y),size(z)) :real, intent(in), optional
: 各微分項の係数
d(size(x),size(y),size(z)) :real, intent(in), optional
: 各微分項の係数
e(size(x),size(y),size(z)) :real, intent(in), optional
: 各微分項の係数
f(size(x),size(y),size(z)) :real, intent(in), optional
: 各微分項の係数
g(size(x),size(y),size(z)) :real, intent(in), optional
: 各微分項の係数
undef :real, intent(in), optional
: 未定義値
inner_bound(size(x),size(y),size(z)) :integer, intent(in), optional
: 内部領域の境界. 値に応じてその格子点で境界値計算 1 = 固定端境界, 10 = 境界の内側. 2 = y 方向自由端境界 (フラックスは上向き) -2 = y 方向自由端境界 (フラックスは下向き) 4 = x 方向自由端境界 (フラックスは右向き) -4 = x 方向自由端境界 (フラックスは左向き) 3 = 周期境界 8 = |_, ~| で両方とも自由境界条件 -8 = |~, _| で両方とも自由境界条件 この引数が与えられなければ全領域を計算する. 境界の内側格子点 (10) は反復計算を行わず, undef で設定された値もしくはゼロが入る. このときの境界値は bound_opt の値が用いられる.
init_flag :logical, intent(in), optional
: psi の値をゼロで初期化するか. .true. = 初期化する. .false. = 初期化しない. デフォルトでは初期化する.
accel :real, intent(in), optional
: SOR の加速係数 (0 < accel < 2) デフォルト = 1
ln :integer, intent(in), optional
: 反復回数 この値が与えられるとき, eps の値に無関係に ln 回ループさせる.

[Source]

subroutine Ellip_Jacobi_3d(x, y, z, rho, eps, boundary, psi, bound_opt, xa, ya, za, a, b, c, d, e, f, g, undef, inner_bound, init_flag, accel, ln )
use omp_lib
! openmp によるスレッド並列が可能.
! ガウス・ザイデル法ではアルゴリズムの点から並列化が困難と思われたので,
! 並列計算によるポアソン方程式の求積が必要となるなら,
! ヤコビ法のものを使用されたい.
! 各オプション配列は, ポアソン系の各微分項の値. デフォルトはゼロで設定される.
! $$xa\dfrac{\partial ^2\psi}{\partial x^2} +ya\dfrac{\partial ^2\psi}{\partial y^2} +za\dfrac{\partial ^2\psi}{\partial z^2} +a\dfrac{\partial ^2\psi}{\partial x\partial y} +b\dfrac{\partial ^2\psi}{\partial y\partial z} +c\dfrac{\partial ^2\psi}{\partial z\partial x} +d\dfrac{\partial \psi}{\partial x} +e\dfrac{\partial \psi}{\partial y} +f\dfrac{\partial \psi}{\partial z} +g\psi =\rho $$
! の各係数に対応している.
  implicit none
  real, intent(in) :: x(:)  ! 領域の x 座標
  real, intent(in) :: y(:)  ! 領域の y 座標
  real, intent(in) :: z(:)  ! 領域の z 座標
  real, intent(in) :: rho(size(x),size(y),size(z))  ! ポアソン方程式の強制項
                   ! rho =0 でラプラス方程式も求積可能
  real, intent(in) :: eps  ! 収束条件
  character(6), intent(in) :: boundary  ! 境界条件
                ! 4 文字で各辺の境界条件を与える.
                ! 1 文字目 : x 下端, 2 文字目 : y 左端, 3 文字目 : x 上端,
                ! 4 文字目 : y 右端
                ! boundary は 1 : 固定端境界, 2 : 自由端境界, 3 : 周期境界
  real, intent(in), optional :: bound_opt(size(x),size(y),size(z))
                             ! 境界での強制
                             ! ノイマン境界の場合 : フラックス値
  real, intent(in), optional :: xa(size(x),size(y),size(z))  ! 各微分項の係数
  real, intent(in), optional :: ya(size(x),size(y),size(z))  ! 各微分項の係数
  real, intent(in), optional :: za(size(x),size(y),size(z))  ! 各微分項の係数
  real, intent(in), optional :: a(size(x),size(y),size(z))  ! 各微分項の係数
  real, intent(in), optional :: b(size(x),size(y),size(z))  ! 各微分項の係数
  real, intent(in), optional :: c(size(x),size(y),size(z))  ! 各微分項の係数
  real, intent(in), optional :: d(size(x),size(y),size(z))  ! 各微分項の係数
  real, intent(in), optional :: e(size(x),size(y),size(z))  ! 各微分項の係数
  real, intent(in), optional :: f(size(x),size(y),size(z))  ! 各微分項の係数
  real, intent(in), optional :: g(size(x),size(y),size(z))  ! 各微分項の係数
  real, intent(inout) :: psi(size(x),size(y),size(z))  ! ポアソン方程式の解
  real, intent(in), optional :: undef  ! 未定義値
  integer, intent(in), optional :: inner_bound(size(x),size(y),size(z))
                             ! 内部領域の境界. 値に応じてその格子点で境界値計算
                             ! 1 = 固定端境界, 10 = 境界の内側.
                             ! 2 = y 方向自由端境界 (フラックスは上向き)
                             ! -2 = y 方向自由端境界 (フラックスは下向き)
                             ! 4 = x 方向自由端境界 (フラックスは右向き)
                             ! -4 = x 方向自由端境界 (フラックスは左向き)
                             ! 3 = 周期境界
                             ! 8 = |_, ~| で両方とも自由境界条件
                             ! -8 = |~, _| で両方とも自由境界条件
                             ! この引数が与えられなければ全領域を計算する.
                             ! 境界の内側格子点 (10) は反復計算を行わず,
                             ! undef で設定された値もしくはゼロが入る.
                             ! このときの境界値は bound_opt の値が用いられる.
  logical, intent(in), optional :: init_flag  ! psi の値をゼロで初期化するか.
                             ! .true. = 初期化する. .false. = 初期化しない.
                             ! デフォルトでは初期化する.
  real, intent(in), optional :: accel  ! SOR の加速係数 (0 < accel < 2)
                             ! デフォルト = 1
  integer, intent(in), optional :: ln  ! 反復回数
                             ! この値が与えられるとき, eps の値に無関係に
                             ! ln 回ループさせる.

  integer :: i, j, k, ix, jy, kz, nl, counter
  integer :: nx  ! x 方向の配列要素
  integer :: ny  ! y 方向の配列要素
  integer :: nz  ! z 方向の配列要素
  integer :: signa, signb, signc  ! 各係数を計算するかどうか
  integer, dimension(size(x),size(y),size(z)) :: ib

  real :: defun
  real :: err, err_max, accc
  real :: bnd(size(x),size(y),size(z))
  real :: dx(size(x)), dy(size(y)), dx2(size(x)), dy2(size(y))
  real :: dz(size(z)), dz2(size(z))
  real, dimension(size(x),size(y)) :: dxdy
  real, dimension(size(y),size(z)) :: dydz
  real, dimension(size(x),size(z)) :: dxdz
  real, dimension(size(x),size(y),size(z)) :: xt, yt, zt, at, bt, ct, dt, et, ft, gt
  real, dimension(size(x),size(y),size(z)) :: xdp, xdm, yep, yem, zfp, zfm, xyz, divi
  real, dimension(size(x),size(y),size(z)) :: tmp, tmp_b

  character(6) :: bound
  logical :: sor_flag
  logical, dimension(size(x),size(y),size(z)) :: inner_flag

  bound(1:6)=boundary(1:6)

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

!-- 応答関数の初期化

  if(present(init_flag))then
     if(init_flag.eqv..true.)then
        psi = 0.0
     end if
  else
     psi = 0.0
  end if

!-- 内部境界の判別フラグの設定

  if(present(inner_bound))then
     call set_bound_3d( bound, ib, inner_flag, inner_bound )
  else
     call set_bound_3d( bound, ib, inner_flag )
  end if

!-- 領域・内部境界における境界値の設定

  if(present(bound_opt))then
     call setval_bound_3d( ib, bnd, psi, bound_opt )
  else
     call setval_bound_3d( ib, bnd, psi )
  end if

!-- 未定義値の設定

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

!-- 係数の代入
!-- xa, ya, za については, 値が入れられていなければ, 全配列に 1 を代入する.
  if(present(xa))then
     call set_coe_3d( xt, ext=xa )
  else
     call set_coe_3d( xt, def=1.0 )
  end if

  if(present(ya))then
     call set_coe_3d( yt, ext=ya )
  else
     call set_coe_3d( yt, def=1.0 )
  end if

  if(present(za))then
     call set_coe_3d( zt, ext=za )
  else
     call set_coe_3d( zt, def=1.0 )
  end if

  if(present(a))then
     call set_coe_3d( at, ext=a )
     signa=1
  else
     call set_coe_3d( at, def=0.0 )
     signa=0
  end if

  if(present(b))then
     call set_coe_3d( bt, ext=b )
     signb=1
  else
     call set_coe_3d( bt, def=0.0 )
     signb=0
  end if

  if(present(c))then
     signc=1
     call set_coe_3d( ct, ext=c )
  else
     call set_coe_3d( ct, def=0.0 )
     signc=0
  end if

  if(present(d))then
     call set_coe_3d( dt, ext=d )
  else
     call set_coe_3d( dt, def=0.0 )
  end if

  if(present(e))then
     call set_coe_3d( et, ext=e )
  else
     call set_coe_3d( et, def=0.0 )
  end if

  if(present(f))then
     call set_coe_3d( ft, ext=f )
  else
     call set_coe_3d( ft, def=0.0 )
  end if

  if(present(g))then
     call set_coe_3d( gt, ext=g )
  else
     call set_coe_3d( gt, def=0.0 )
  end if

!-- 最高階数における係数チェック. (係数がゼロでないか調べる.)

  call check_coe_3d( xt, 0.0 )
  call check_coe_3d( yt, 0.0 )
  call check_coe_3d( zt, 0.0 )

!-- 加速係数の判定

  if(present(accel))then
     accc=accel
     sor_flag=.true.
  else
     sor_flag=.false.
  end if

!-- ln の処理
  if(present(ln))then
     nl=ln
  else
     nl=0
  end if

!-- 以下で先に格子間隔等の 1 回計算でよいものを求めておく.
!-- これらは 1 方向のみで変化すればよい.
!-- 格子点間隔の計算
  do i=2,nx-1
     dx(i)=(x(i+1)-x(i-1))*0.5
     dx2(i)=dx(i)**2
  end do
  do j=2,ny-1
     dy(j)=(y(j+1)-y(j-1))*0.5
     dy2(j)=dy(j)**2
  end do
  do k=2,nz-1
     dz(k)=(z(k+1)-z(k-1))*0.5
     dz2(k)=dz(k)**2
  end do

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

  do j=1,ny
     do i=1,nx
        dxdy(i,j)=dx(i)*dy(j)
     end do
  end do
  do k=1,nz
     do j=1,ny
        dydz(j,k)=dy(j)*dz(k)
     end do
  end do
  do k=1,nz
     do i=1,nx
        dxdz(i,k)=dx(i)*dz(k)
     end do
  end do

!-- ポアソン係数の計算
!-- 実際にポアソン係数が使用されるのは, 境界を除いた 1 回り内側なので,
!-- 計算量削減のため, ループをそのようにしておく.

  xyz=0.0
  xdp=0.0
  xdm=0.0
  yep=0.0
  yem=0.0
  zfp=0.0
  zfm=0.0

!-- 最高次数係数 ac の計算
!$omp parallel default(shared)
!$omp do schedule(dynamic) private(i,j,k)

  do k=2,nz-1
     do j=2,ny-1
        do i=2,nx-1
           xyz(i,j,k)=1.0/(2.0*(xt(i,j,k)/dx2(i)+yt(i,j,k)/dy2(j)+zt(i,j,k)/dz2(k))-gt(i,j,k))
           xdp(i,j,k)=(xt(i,j,k)/(dx2(i))+0.5*dt(i,j,k)/dx(i))*xyz(i,j,k)
           xdm(i,j,k)=(xt(i,j,k)/(dx2(i))-0.5*dt(i,j,k)/dx(i))*xyz(i,j,k)
           yep(i,j,k)=(yt(i,j,k)/(dy2(j))+0.5*et(i,j,k)/dy(j))*xyz(i,j,k)
           yem(i,j,k)=(yt(i,j,k)/(dy2(j))-0.5*et(i,j,k)/dy(j))*xyz(i,j,k)
           zfp(i,j,k)=(zt(i,j,k)/(dz2(k))+0.5*ft(i,j,k)/dz(k))*xyz(i,j,k)
           zfm(i,j,k)=(zt(i,j,k)/(dz2(k))-0.5*ft(i,j,k)/dz(k))*xyz(i,j,k)
           if(signa==1)then
              at(i,j,k)=0.25*at(i,j,k)/(dxdy(i,j))*xyz(i,j,k)
           else
              at(i,j,k)=0.0
           end if
           if(signb==1)then
              bt(i,j,k)=0.25*bt(i,j,k)/(dydz(j,k))*xyz(i,j,k)
           else
              bt(i,j,k)=0.0
           end if
           if(signc==1)then
              ct(i,j,k)=0.25*ct(i,j,k)/(dxdz(i,k))*xyz(i,j,k)
           else
              ct(i,j,k)=0.0
           end if
        end do
     end do
  end do

!$omp end do
!$omp end parallel

  err_max=eps  ! while に入るための便宜的措置
  counter=0

!-- 実際のソルバ ---
  do while(err_max>=eps)
     err_max=0.0

!$omp parallel default(shared)
!$omp do schedule(dynamic) private(i,j,k,ix,jy,kz)
     do k=1,nz
        do j=1,ny
           do i=1,nx

!-- 以降, 反復計算に必要な点の周囲 8 点についてそれぞれ
!-- inner_flag のチェックと同時に適切な値をそれぞれ逐次計算する,
!-- 下で, 各方向の格子計算を行っているが, ib=8,-8 は 4 隅格子のときにしか
!-- case select しない. なぜなら, 真横や上下に隅格子がくることはありえない.
              if(inner_flag(i,j,k).eqv..false.)then  ! .false. なら領域計算開始

                 tmp(i,j,k)=-rho(i,j,k)*xyz(i,j,k)
                 tmp(i,j,k)=tmp(i,j,k)+xdp(i,j,k)*psi(i+1,j,k) +xdm(i,j,k)*psi(i-1,j,k) +yep(i,j,k)*psi(i,j+1,k) +yem(i,j,k)*psi(i,j-1,k) +zfp(i,j,k)*psi(i,j,k+1) +zfm(i,j,k)*psi(i,j,k-1)


                 if(signa==1)then
                    tmp_b(i,j,k)=at(i,j,k)*(psi(i+1,j+1,k)+psi(i-1,j-1,k) -psi(i+1,j-1,k)-psi(i-1,j+1,k))
                 else
                    tmp_b(i,j,k)=0.0
                 end if
                 if(signb==1)then  ! そもそも bt = 0 なら計算しない.
                    tmp_b(i,j,k)=tmp_b(i,j,k) +bt(i,j,k)*(psi(i,j+1,k+1)+psi(i,j-1,k-1) -psi(i,j-1,k+1)-psi(i,j+1,k-1))
                 end if
                 if(signc==1)then  ! そもそも bt = 0 なら計算しない.
                    tmp_b(i,j,k)=tmp_b(i,j,k) +ct(i,j,k)*(psi(i+1,j,k+1)+psi(i-1,j,k-1) -psi(i-1,j,k+1)-psi(i+1,j,k-1))
                 end if
              
                 tmp(i,j,k)=tmp(i,j,k)+tmp_b(i,j,k)

              else   ! .true. なら境界計算開始.

                 select case (ib(i,j,k))
                 case (1)
                    tmp(i,j,k)=bnd(i,j,k)

                 case (2)  ! x 方向にフラックス一定, 上側が参照値
                    tmp(i,j,k)=psi(i+1,j,k)-bnd(i,j,k)*dx(i)

                 case (-2)  ! x 方向にフラックス一定, 下側が参照値
                    tmp(i,j,k)=psi(i-1,j,k)+bnd(i,j,k)*dx(i)

                 case (4)  ! y 方向にフラックス一定, 右側が参照値
                    tmp(i,j,k)=psi(i,j+1,k)-bnd(i,j,k)*dy(j)

                 case (-4)  ! y 方向にフラックス一定, 左側が参照値
                    tmp(i,j,k)=psi(i,j-1,k)+bnd(i,j,k)*dy(j)

                 case (6)  ! y 方向にフラックス一定, 右側が参照値
                    tmp(i,j,k)=psi(i,j,k+1)-bnd(i,j,k)*dz(k)

                 case (-6)  ! y 方向にフラックス一定, 左側が参照値
                    tmp(i,j,k)=psi(i,j,k-1)+bnd(i,j,k)*dz(k)

                 case (3)  ! 12 辺, もしくは 8 点で周期境界を判断
                    if(i==1)then
                       ix=nx-1
                    else if(i==nx)then
                       ix=2
                    else
                       ix=i
                    end if
                    if(j==1)then
                       jy=ny-1
                    else if(j==ny)then
                       jy=2
                    else
                       jy=j
                    end if
                    if(k==1)then
                       kz=nz-1
                    else if(k==nz)then
                       kz=2
                    else
                       kz=k
                    end if
                    tmp(i,j,k)=psi(ix,jy,kz)

                 case (8)  ! 両方フラックス一定で z 面の x, y 右上か左下角.
                    if(i==1.and.j==1)then  ! -- 評価 1
                       tmp(i,j,k)=psi(i+1,j+1,k) -0.5*(bnd(i+1,j,k)*dy(j)+bnd(i,j+1,k)*dx(i))

                    else if(i==nx.and.j==ny)then  ! -- 評価 2
                       tmp(i,j,k)=psi(i-1,j-1,k) +0.5*(bnd(i-1,j,k)*dy(j)+bnd(i,j-1,k)*dx(i))

                    else if(ib(i-1,j,k)==10.and.ib(i,j-1,k)==10)then
                       ! -- 評価 1 と同じ
                       tmp(i,j,k)=psi(i+1,j+1,k) -0.5*(bnd(i+1,j,k)*dy(j)+bnd(i,j+1,k)*dx(i))

                    else if(ib(i+1,j,k)==10.and.ib(i,j+1,k)==10)then
                       ! -- 評価 2 と同じ
                       tmp(i,j,k)=psi(i-1,j-1,k) +0.5*(bnd(i-1,j,k)*dy(j)+bnd(i,j-1,k)*dx(i))

                    end if

                 case (-8)  ! 両方フラックス一定で z 面の x, y 右下か左上角.
                    if(i==1.and.j==ny)then  ! -- 評価 1
                       tmp(i,j,k)=psi(i+1,j-1,k) +0.5*(bnd(i+1,j,k)*dy(j)-bnd(i,j-1,k)*dx(i))

                    else if(i==nx.and.j==1)then  ! -- 評価 2
                       tmp(i,j,k)=psi(i-1,j+1,k) +0.5*(-bnd(i-1,j,k)*dy(j)+bnd(i,j+1,k)*dx(i))

                    else if(ib(i-1,j,k)==10.and.ib(i,j+1,k)==10)then
                       ! -- 評価 1 と同じ
                       tmp(i,j,k)=psi(i+1,j-1,k) +0.5*(bnd(i+1,j,k)*dy(j)-bnd(i,j-1,k)*dx(i))

                    else if(ib(i+1,j,k)==10.and.ib(i,j-1,k)==10)then
                       ! -- 評価 2 と同じ
                       tmp(i,j,k)=psi(i-1,j+1,k) +0.5*(-bnd(i-1,j,k)*dy(j)+bnd(i,j+1,k)*dx(i))
                    end if

                 case (12)  ! 両方フラックス一定で y 面の x, z 右上か左下角.
                    if(i==1.and.k==1)then  ! -- 評価 1
                       tmp(i,j,k)=psi(i+1,j,k+1) -0.5*(bnd(i+1,j,k)*dz(k)+bnd(i,j,k+1)*dx(i))

                    else if(i==nx.and.k==nz)then  ! -- 評価 2
                       tmp(i,j,k)=psi(i-1,j,k-1) +0.5*(bnd(i-1,j,k)*dz(k)+bnd(i,j,k-1)*dx(i))

                    else if(ib(i-1,j,k)==10.and.ib(i,j,k-1)==10)then
                       ! -- 評価 1 と同じ
                       tmp(i,j,k)=psi(i+1,j,k+1) -0.5*(bnd(i+1,j,k)*dz(k)+bnd(i,j,k+1)*dx(i))

                    else if(ib(i+1,j,k)==10.and.ib(i,j,k+1)==10)then
                       ! -- 評価 2 と同じ
                       tmp(i,j,k)=psi(i-1,j,k-1) +0.5*(bnd(i-1,j,k)*dz(k)+bnd(i,j,k-1)*dx(i))

                    end if

                 case (-12)  ! 両方フラックス一定で y 面の x, z 右下か左上角.
                    if(i==1.and.k==nz)then  ! -- 評価 1
                       tmp(i,j,k)=psi(i+1,j,k-1) +0.5*(bnd(i+1,j,k)*dz(k)-bnd(i,j,k-1)*dx(i))

                    else if(i==nx.and.k==1)then  ! -- 評価 2
                       tmp(i,j,k)=psi(i-1,j,k+1) +0.5*(-bnd(i-1,j,k)*dz(k)+bnd(i,j,k+1)*dx(i))

                    else if(ib(i-1,j,k)==10.and.ib(i,j,k+1)==10)then
                       ! -- 評価 1 と同じ
                       tmp(i,j,k)=psi(i+1,j,k-1) +0.5*(bnd(i+1,j,k)*dz(k)-bnd(i,j,k-1)*dx(i))

                    else if(ib(i+1,j,k)==10.and.ib(i,j,k-1)==10)then
                       ! -- 評価 2 と同じ
                       tmp(i,j,k)=psi(i-1,j,k+1) +0.5*(-bnd(i-1,j,k)*dz(k)+bnd(i,j,k+1)*dx(i))

                    end if

                 case (24)  ! 両方フラックス一定で x 面の y, z 右上か左下角.
                    if(j==1.and.k==1)then  ! -- 評価 1
                       tmp(i,j,k)=psi(i,j+1,k+1) -0.5*(bnd(i,j+1,k)*dz(k)+bnd(i,j,k+1)*dy(j))

                    else if(j==ny.and.k==nz)then  ! -- 評価 2
                       tmp(i,j,k)=psi(i,j-1,k-1) +0.5*(bnd(i,j-1,k)*dz(k)+bnd(i,j,k-1)*dy(j))

                    else if(ib(i,j-1,k)==10.and.ib(i,j,k-1)==10)then
                       ! -- 評価 1 と同じ
                       tmp(i,j,k)=psi(i,j+1,k+1) -0.5*(bnd(i,j+1,k)*dz(k)+bnd(i,j,k+1)*dy(j))

                    else if(ib(i,j+1,k)==10.and.ib(i,j,k+1)==10)then
                       ! -- 評価 2 と同じ
                       tmp(i,j,k)=psi(i,j-1,k-1) +0.5*(bnd(i,j-1,k)*dz(k)+bnd(i,j,k-1)*dy(j))

                    end if

                 case (-24)  ! 両方フラックス一定で x 面の y, z 右下か左上角.
                    if(j==1.and.k==nz)then  ! -- 評価 1
                       tmp(i,j,k)=psi(i,j+1,k-1) +0.5*(bnd(i,j+1,k)*dz(k)-bnd(i,j,k-1)*dy(j))

                    else if(j==ny.and.k==1)then  ! -- 評価 2
                       tmp(i,j,k)=psi(i,j-1,k+1) +0.5*(-bnd(i,j-1,k)*dz(k)+bnd(i,j,k+1)*dy(j))

                    else if(ib(i,j-1,k)==10.and.ib(i,j,k+1)==10)then
                       ! -- 評価 1 と同じ
                       tmp(i,j,k)=psi(i,j+1,k-1) +0.5*(bnd(i,j+1,k)*dz(k)-bnd(i,j,k-1)*dy(j))

                    else if(ib(i,j+1,k)==10.and.ib(i,j,k-1)==10)then
                       ! -- 評価 2 と同じ
                       tmp(i,j,k)=psi(i,j-1,k+1) +0.5*(-bnd(i,j-1,k)*dz(k)+bnd(i,j,k+1)*dy(j))

                    end if

                 !-- 以降, 隅領域なので, 個別に設定.
                 case (11)  ! 両方フラックス一定で (1,1,1) 点, もしくは内部領域.
                    tmp(i,j,k)=psi(i+1,j+1,k+1) -(bnd(i,j+1,k+1)*dx(i) +bnd(i+1,j,k+1)*dy(j) +bnd(i+1,j+1,k)*dz(k))/3.0

                 case (13)  ! 両方フラックス一定で (nx,1,1) 点, もしくは内部領域.
                    tmp(i,j,k)=psi(i-1,j+1,k+1) -(-bnd(i,j+1,k+1)*dx(i) +bnd(i-1,j,k+1)*dy(j) +bnd(i-1,j+1,k)*dz(k))/3.0

                 case (17)  ! 両方フラックス一定で (1,ny,1) 点, もしくは内部領域.
                    tmp(i,j,k)=psi(i+1,j-1,k+1) -(bnd(i,j-1,k+1)*dx(i) -bnd(i+1,j,k+1)*dy(j) +bnd(i+1,j-1,k)*dz(k))/3.0

                 case (19)  ! 両方フラックス一定で (nx,ny,1) 点, もしくは内部領域.
                    tmp(i,j,k)=psi(i-1,j-1,k+1) -(-bnd(i,j-1,k+1)*dx(i) -bnd(i-1,j,k+1)*dy(j) +bnd(i-1,j-1,k)*dz(k))/3.0

                 case (23)  ! 両方フラックス一定で (1,1,nz) 点, もしくは内部領域.
                    tmp(i,j,k)=psi(i+1,j+1,k-1) -(bnd(i,j+1,k-1)*dx(i) +bnd(i+1,j,k-1)*dy(j) -bnd(i+1,j+1,k)*dz(k))/3.0

                 case (29)  ! 両方フラックス一定で (nx,1,nz) 点, もしくは内部領域.
                    tmp(i,j,k)=psi(i-1,j+1,k-1) -(-bnd(i,j+1,k-1)*dx(i) +bnd(i-1,j,k-1)*dy(j) -bnd(i-1,j+1,k)*dz(k))/3.0

                 case (31)  ! 両方フラックス一定で (1,ny,nz) 点, もしくは内部領域.
                    tmp(i,j,k)=psi(i+1,j-1,k-1) -(bnd(i,j-1,k-1)*dx(i) -bnd(i+1,j,k-1)*dy(j) -bnd(i+1,j-1,k)*dz(k))/3.0

                 case (37)  ! 両方フラックス一定で (nx,ny,nz) 点, もしくは内部領域.
                    tmp(i,j,k)=psi(i-1,j-1,k-1) -(-bnd(i,j-1,k-1)*dx(i) -bnd(i-1,j,k-1)*dy(j) -bnd(i-1,j-1,k)*dz(k))/3.0

                 case (-11)  ! 両方フラックス一定で (i+1,j+1,k+1) に 10 が設定.
                    tmp(i,j,k)=(psi(i-1,j+1,k+1)+psi(i+1,j-1,k+1)+psi(i+1,j+1,k-1) +bnd(i,j+1,k+1)*dx(i)+bnd(i+1,j,k+1)*dy(j) +bnd(i+1,j+1,k)*dz(k))/3.0

                 case (-13)  ! 両方フラックス一定で (i-1,j+1,k+1) に 10 が設定.
                    tmp(i,j,k)=(psi(i+1,j+1,k+1)+psi(i-1,j-1,k+1)+psi(i-1,j+1,k-1) -bnd(i,j+1,k+1)*dx(i)+bnd(i-1,j,k+1)*dy(j) +bnd(i-1,j+1,k)*dz(k))/3.0

                 case (-17)  ! 両方フラックス一定で (i+1,j-1,k+1) に 10 が設定.
                    tmp(i,j,k)=(psi(i-1,j-1,k+1)+psi(i+1,j+1,k+1)+psi(i+1,j-1,k-1) +bnd(i,j-1,k+1)*dx(i)-bnd(i+1,j,k+1)*dy(j) +bnd(i+1,j-1,k)*dz(k))/3.0

                 case (-19)  ! 両方フラックス一定で (i-1,j-1,k+1) に 10 が設定.
                    tmp(i,j,k)=(psi(i+1,j-1,k+1)+psi(i-1,j+1,k+1)+psi(i-1,j-1,k-1) -bnd(i,j-1,k+1)*dx(i)-bnd(i-1,j,k+1)*dy(j) +bnd(i-1,j-1,k)*dz(k))/3.0

                 case (-23)  ! 両方フラックス一定で (i+1,j+1,k-1) に 10 が設定.
                    tmp(i,j,k)=(psi(i-1,j+1,k-1)+psi(i+1,j-1,k-1)+psi(i+1,j+1,k+1) +bnd(i,j+1,k-1)*dx(i)+bnd(i+1,j,k-1)*dy(j) -bnd(i+1,j+1,k)*dz(k))/3.0

                 case (-29)  ! 両方フラックス一定で (i-1,j+1,k-1) に 10 が設定.
                    tmp(i,j,k)=(psi(i+1,j+1,k-1)+psi(i-1,j-1,k-1)+psi(i-1,j+1,k+1) -bnd(i,j+1,k-1)*dx(i)+bnd(i-1,j,k-1)*dy(j) -bnd(i-1,j+1,k)*dz(k))/3.0

                 case (-31)  ! 両方フラックス一定で (i+1,j-1,k-1) に 10 が設定.
                    tmp(i,j,k)=(psi(i-1,j-1,k-1)+psi(i+1,j+1,k-1)+psi(i+1,j-1,k+1) +bnd(i,j-1,k-1)*dx(i)-bnd(i+1,j,k-1)*dy(j) -bnd(i+1,j-1,k)*dz(k))/3.0

                 case (-37)  ! 両方フラックス一定で (i-1,j-1,k-1) に 10 が設定.
                    tmp(i,j,k)=(psi(i+1,j-1,k-1)+psi(i-1,j+1,k-1)+psi(i-1,j-1,k+1) -bnd(i,j-1,k-1)*dx(i)-bnd(i-1,j,k-1)*dy(j) -bnd(i-1,j-1,k)*dz(k))/3.0

                 end select

              end if

              if(sor_flag.eqv..true.)then
                 tmp(i,j,k)=(1.0-accc)*psi(i,j,k)+tmp(i,j,k)*accc
              end if

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

!-- 最大誤差の更新
     do k=2,nz-1
        do j=2,ny-1
           do i=2,nx-1
              err=abs(tmp(i,j,k)-psi(i,j,k))
              if(err_max<=err)then
                 err_max=err
              end if
           end do
        end do
     end do

!-- 一斉更新

     do k=1,nz
        do j=1,ny
           do i=1,nx
              psi(i,j,k)=tmp(i,j,k)
           end do
        end do
     end do

     if(nl/=0)then
        counter=counter+1
        if(counter==nl)then
           exit
        else
           err_max=eps
        end if
     end if

  end do

!-- 境界の設定

  call calculate_bound_3d( ib, dx, dy, dz, bnd, psi )

!-- 未定義領域には undef を代入する.

  do k=1,nz
     do j=1,ny
        do i=1,nx
           if(ib(i,j,k)==10)then
              psi(i,j,k)=defun
           end if
        end do
     end do
  end do

end subroutine Ellip_Jacobi_3d
Subroutine :
coe(:,:) :real, intent(inout)
: 代入される配列
aval :real, intent(in)
: チェックされる値

2 次元配列に aval が代入されていないかをチェックする.

[Source]

subroutine check_coe( coe, aval )
! 2 次元配列に aval が代入されていないかをチェックする.
  implicit none
  real, intent(inout) :: coe(:,:)  ! 代入される配列
  real, intent(in) :: aval         ! チェックされる値
  integer :: i, j, nx, ny

  nx=size(coe,1)
  ny=size(coe,2)

  do j=1,ny
     do i=1,nx
        if(coe(i,j)==aval)then
           write(*,*) "### ERROR (Ellip_Slv module) ###"
           write(*,*) "Detect a certain value", aval
           write(*,*) "STOP."
           stop
        end if
     end do
  end do

end subroutine check_coe
Subroutine :
coe(:,:,:) :real, intent(inout)
: 代入される配列
aval :real, intent(in)
: 検出される値

3 次元配列に aval で指定された値が入っていないかを検出する.

[Source]

subroutine check_coe_3d( coe, aval )
! 3 次元配列に aval で指定された値が入っていないかを検出する.
  implicit none
  real, intent(inout) :: coe(:,:,:)  ! 代入される配列
  real, intent(in) :: aval           ! 検出される値
  integer :: i, j, k, nx, ny, nz

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

  do k=1,nz
     do j=1,ny
        do i=1,nx
           if(coe(i,j,k)==aval)then
              write(*,*) "### ERROR (Ellip_Slv module) ###"
              write(*,*) "Detect a certain value", aval
              write(*,*) "STOP."
              stop
           end if
        end do
     end do
  end do

end subroutine check_coe_3d