Class Alge_Solv
In: alge_solv.f90

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

Methods

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
: 各微分項の係数 ノイマン境界の場合 : フラックス値

ガウス=ザイデル法によるポアソン方程式の求積 各オプション配列は, ポアソン系の各微分項の値. デフォルトはゼロで設定される. $$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} =rho $$ の各係数に対応している.

[Source]

subroutine Poisson_GauSei(x, y, rho, eps, boundary, psi, bound_opt, a, b, c, d, e)
! ガウス=ザイデル法によるポアソン方程式の求積
! 各オプション配列は, ポアソン系の各微分項の値. デフォルトはゼロで設定される.
! $$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} =\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(inout) :: psi(size(x),size(y))  ! ポアソン方程式の解
  integer :: i, j, k, l, m, n
  integer :: nx  ! x 方向の配列要素
  integer :: ny  ! y 方向の配列要素
  integer :: signb, signd, signe  ! 各係数を計算するかどうか
  real :: tmp, err, err_max
  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
  real, dimension(size(x),size(y)) :: adp, adm, cep, cem, ac
  character(4) :: bound
  real :: tmp_b, tmp_adm, tmp_adp, tmp_cem, tmp_cep

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

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

  psi = 0.0

!-- 周期境界の設定確認.
!-- 周期境界なので, 両端とも 3 が設定されていないといけない.
  if(bound(1:1)=='3')then
     if(bound(3:3)/='3')then
        write(*,*) "### ERROR ###"
        write(*,*) "if bound = 3, bound(1:1)==bound(3:3). STOP."
        stop
     end if
  end if

  if(bound(3:3)=='3')then
     if(bound(1:1)/='3')then
        write(*,*) "### ERROR ###"
        write(*,*) "if bound = 3, bound(1:1)==bound(3:3). STOP."
        stop
     end if
  end if

  if(bound(2:2)=='3')then
     if(bound(4:4)/='3')then
        write(*,*) "### ERROR ###"
        write(*,*) "if bound = 3, bound(2:2)==bound(4:4). STOP."
        stop
     end if
  end if

  if(bound(4:4)=='3')then
     if(bound(2:2)/='3')then
        write(*,*) "### ERROR ###"
        write(*,*) "if bound = 3, bound(2:2)==bound(4:4). STOP."
        stop
     end if
  end if

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

  if(present(c))then
     do j=1,ny
        do i=1,nx
           ct(i,j)=c(i,j)
        end do
     end do
  else
     do j=1,ny
        do i=1,nx
           ct(i,j)=1.0
        end do
     end do
  end if

  if(present(b))then
     do j=1,ny
        do i=1,nx
           bt(i,j)=b(i,j)
        end do
     end do
     signb=1
  else
     signb=0
  end if

  if(present(d))then
     do j=1,ny
        do i=1,nx
           dt(i,j)=d(i,j)
        end do
     end do
     signd=1
  else
     signd=0
  end if

  if(present(e))then
     do j=1,ny
        do i=1,nx
           et(i,j)=e(i,j)
        end do
     end do
     signe=1
  else
     signe=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

!-- ポアソン係数の計算

  if(signd==0)then  ! 付加項 d がついていないとき
     do j=2,ny-1
        do i=2,nx-1
           adp(i,j)=at(i,j)/(dx2(i))
           adm(i,j)=at(i,j)/(dx2(i))
        end do
     end do
  else
     do j=2,ny-1
        do i=2,nx-1
           adp(i,j)=at(i,j)/(dx2(i))+0.5*dt(i,j)/dx(i)
           adm(i,j)=at(i,j)/(dx2(i))-0.5*dt(i,j)/dx(i)
        end do
     end do
  end if

  if(signe==0)then  ! 付加項 e がついていないとき
     do j=2,ny-1
        do i=2,nx-1
           cep(i,j)=ct(i,j)/(dy2(j))
           cem(i,j)=ct(i,j)/(dy2(j))
        end do
     end do
  else
     do j=2,ny-1
        do i=2,nx-1
           cep(i,j)=ct(i,j)/(dy2(j))+0.5*et(i,j)/dy(j)
           cem(i,j)=ct(i,j)/(dy2(j))-0.5*et(i,j)/dy(j)
        end do
     end do
  end if

!-- 最高次数係数 ac の計算 (境界条件によって評価式が変わる.)
!-- 実際にポアソン係数が使用されるのは, 境界を除いた 1 回り内側なので,
!-- 計算量削減のため, ループをそのようにしておく.
  do j=2,ny-1  ! 係数 ac については, 境界の内側で境界条件によって式が異なる.
     do i=2,nx-1  ! ここでは, 固定端条件での値を代入しており, それ以外の場合は以下の処理で上書きする.
        ac(i,j)=2.0*(at(i,j)/dx2(i)+ct(i,j)/dy2(j))
     end do
  end do

!-- 以下, 境界の内側での係数 ac の計算
!-- 本計算では, 固定境界があるなら, 隅領域はすべてその値で固定するようにしており,
!-- 固定境界値は反復計算の前に与えられる値を用いるので, 境界条件 1 では評価式は変わらない.
!-- また, 周期境界条件では, 両端が周期境界でなければならないから, bound の組み合わせは
!-- 制限される.
!-- 実際に係数が変わるのは, 条件 2 のときだけであるので, 2 について評価式を計算する.
!-- (3 は係数ではなく, ポアソンの式だけが変化する.)
!-- 以下では, 境界の 4 隅の内側も計算している.
!-- これで 2 重勘定されるのは, 隅の両側ともが 2 の条件のときだけなので,
!-- この方法を用いた方が場合分けが少なくて済む.
!-- 2 重勘定したときは, その半分を引く.

  if(bound(1:1)=='2')then  ! x 下端がノイマン条件 : psi(i,1)=psi(i,2)-f(i,1)*dy(1)
     do i=2,nx-1
        ac(i,2)=ac(i,2)-cem(i,2)
     end do
  end if

  if(bound(2:2)=='2')then  ! y 左端がノイマン条件 : psi(1,j)=psi(2,j)-f(1,j)*dx(1)
     do j=2,ny-1
        ac(2,j)=ac(2,j)-adm(2,j)
     end do
  end if

  if(bound(3:3)=='2')then  ! x 上端がノイマン条件 : psi(i,ny)=psi(i,ny-1)+f(i,ny)*dy(ny)
     do i=2,nx-1
        ac(i,ny-1)=ac(i,ny-1)-cep(i,ny-1)
     end do
  end if

  if(bound(4:4)=='2')then  ! y 右端がノイマン条件 : psi(nx,j)=psi(nx-1,j)+f(nx,j)*dx(nx)
     do j=2,ny-1
        ac(nx-1,j)=ac(nx-1,j)-adp(nx-1,j)
     end do
  end if

!-- 内側 4 隅での 2 重勘定の解消  bt についてはここでのみ計算されるので, 足し合わせておく.

  if(bound(1:2)=='22')then  ! 左下隅の解消
     ac(2,2)=ac(2,2)+(cem(2,2)+adm(2,2))*0.5-0.25*bt(2,2)/dxdy(2,2)
  end if

  if(bound(2:3)=='22')then  ! 左上隅の解消
     ac(2,ny-1)=ac(2,ny-1)+(cep(2,ny-1)+adm(2,ny-1))*0.5+0.25*bt(2,ny-1)/dxdy(2,ny-1)
  end if

  if(bound(3:4)=='22')then  ! 右上隅の解消
     ac(nx-1,ny-1)=ac(nx-1,ny-1)+(cep(nx-1,ny-1)+adp(nx-1,ny-1))*0.5 -0.25*bt(nx-1,ny-1)/dxdy(nx-1,ny-1)
  end if

  if(bound(1:1)=='2'.and.bound(4:4)=='2')then  ! 右下隅の解消
     ac(nx-1,2)=ac(nx-1,2)+(cem(nx-1,2)+adp(nx-1,2))*0.5+0.25*bt(nx-1,2)/dxdy(nx-1,2)
  end if

!-- 境界値の設定
  if(present(bound_opt))then
     do j=1,ny
        do i=1,nx
           bnd(i,j)=bound_opt(i,j)
        end do
     end do
  else
     do j=1,ny
        do i=1,nx
           bnd(i,j)=0.0
        end do
     end do
  end if

!-- 境界条件の確認 (+ 固定境界なら, その値を代入する.)
  if(bound(1:1)=='1')then
     do i=1,nx
        psi(i,1)=bnd(i,1)
     end do
  end if

  if(bound(2:2)=='1')then
     do j=1,ny
        psi(1,j)=bnd(1,j)
     end do
  end if

  if(bound(3:3)=='1')then
     do i=1,nx
        psi(i,ny)=bnd(i,ny)
     end do
  end if

  if(bound(4:4)=='1')then
     do j=1,ny
        psi(1,j)=bnd(1,j)
     end do
  end if

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

!-- 実際のソルバ ---
  do while(err_max>=eps)
     err_max=0.0
  do j=2,ny-1
     do i=2,nx-1

        if(i/=2.and.i/=nx-1.and.j/=2.and.j/=nx-1)then
           ! 境界の 1 つ内側以外
           tmp=adp(i,j)*psi(i+1,j)+cep(i,j)*psi(i,j+1) +adm(i,j)*psi(i-1,j)+cem(i,j)*psi(i,j-1) -rho(i,j)

           if(signb==1)then
              tmp=tmp+0.25*bt(i,j)*(psi(i+1,j+1)+psi(i-1,j-1) -psi(i+1,j-1)-psi(i-1,j+1))/dxdy(i,j)
           end if

           tmp=tmp/ac(i,j)

        else

!-- 以下で境界の内側計算
!-- コードを見やすくするため, adm, adp, cem, cep の項ごとに case 分けする.

           if(i==2)then  ! 左境界での扱い. ここでは, adm 以外の項は共通計算できる.
              if(j/=2.or.j/=ny-1)then  ! 隅境界以外
                 tmp=adp(i,j)*psi(i+1,j)+cep(i,j)*psi(i,j+1) +cem(i,j)*psi(i,j-1)-rho(i,j)

                 select case (bound(2:2))  ! y 軸左端
                 case ('1')
                    tmp_adm=adm(i,j)*psi(i-1,j)

                    if(signb==1)then
                       tmp_b=0.25*bt(i,j)*(psi(i+1,j+1)+psi(i-1,j-1) -psi(i+1,j-1)-psi(i-1,j+1))/dxdy(i,j)
                    else
                       tmp_b=0.0
                    end if

                 case ('2')  ! u(i-1,j)=u(i,j)-f(i-1,j)*dx(i-1)
                    tmp_adm=-adm(i,j)*bnd(i-1,j)*dx(i-1)
            
                    if(signb==1)then
                       tmp_b=0.25*bt(i,j)*(psi(i+1,j+1) +psi(i,j-1)-bnd(i-1,j-1)*dx(i-1) -psi(i+1,j-1)-psi(i,j+1)+bnd(i-1,j+1)*dx(i-1))/dxdy(i,j)
                    else
                       tmp_b=0.0
                    end if

                 case ('3')  ! u(1,j)=u(nx-1,j)
                    tmp_adm=adm(i,j)*psi(nx-1,j)

                    if(signb==1)then
                       tmp_b=0.25*bt(i,j)*(psi(i+1,j+1)+psi(nx-1,j-1) -psi(i+1,j-1)-psi(nx-1,j+1))/dxdy(i,j)
                    else
                       tmp_b=0.0
                    end if

                 end select

                 tmp=(tmp+tmp_adm+tmp_b)/ac(i,j)

              else

                 if(j==2)then  ! 左下隅  ! adp, cep は境界条件に依存しない.
                    tmp=adp(i,j)*psi(i+1,j)+cep(i,j)*psi(i,j+1) -rho(i,j)

                    select case (bound(1:1))  ! x 軸下端条件 -> cem 項の計算
                    case ('1')
                       tmp_cem=cem(i,j)*psi(i,j-1)

                    case ('2')  ! u(i,j-1)=u(i,j)-f(i,j-1)*dy(j-1)
                       tmp_cem=-cem(i,j)*bnd(i,j-1)*dy(j-1)

                    case ('3')  ! u(i,1)=u(i,ny-1)
                       tmp_cem=cem(i,j)*psi(i,ny-1)

                    end select

                    select case (bound(2:2))
                    case ('1')
                       tmp_adm=adm(i,j)*psi(i-1,j)

                    case ('2')  ! u(i-1,j)=u(i,j)-f(i-1,j)*dx(i-1)
                       tmp_adm=-adm(i,j)*bnd(i-1,j)*dx(i-1)

                    case ('3')  ! u(1,j)=u(nx-1,j)
                       tmp_adm=adm(i,j)*psi(nx-1,j)

                    end select

                    if(signb==1)then
                       select case (bound(1:2))
                       case ('11')  ! 両方固定境界
                          tmp_b=0.25*bt(i,j)*(psi(i+1,j+1)+psi(i-1,j-1) -psi(i+1,j-1)-psi(i-1,j+1))/dxdy(i,j)
                       case ('12')  ! x 固定 y ノイマン
                          ! u(i-1,j)=u(i,j)-f(i-1,j)*dx(i-1), j-1 は固定強制.
                          tmp_b=0.25*bt(i,j)*(psi(i+1,j+1)+psi(i-1,j-1) -psi(i+1,j-1)-psi(i,j+1)+bnd(i-1,j+1)*dx(i-1))/dxdy(i,j)
                       case ('13')  ! x 固定 y 周期  u(1,j)=u(nx-1,j)
                          tmp_b=0.25*bt(i,j)*(psi(i+1,j+1)+psi(i-1,j-1) -psi(i+1,j-1)-psi(nx-1,j+1))/dxdy(i,j)
                       case ('21')  ! x ノイマン y 固定  u(i,j-1)=u(i,j)-f(i,j-1)*dy(j-1)
                          tmp_b=0.25*bt(i,j)*(psi(i+1,j+1)+psi(i-1,j-1) -psi(i+1,j)+bnd(i+1,j-1)*dy(j-1)-psi(i-1,j+1))/dxdy(i,j)
                       case ('22')  ! x, y ノイマン  u(i-1,j)=u(i,j)-f(i-1,j)*dx(i-1)
                          ! u(i-1,j-1)=u(i,j)-0.5*(f(i-1,j)*dx(i-1)+f(i,j-1)*dy(j-1))
                          ! u(i,j-1)=u(i,j)-f(i,j-1)*dy(j-1)
                          tmp_b=0.25*bt(i,j)*(psi(i+1,j+1) -0.5*(bnd(i-1,j)*dx(i-1)+bnd(i,j-1)*dy(j-1)) -psi(i+1,j)+bnd(i+1,j-1)*dy(j-1) -psi(i,j+1)+bnd(i-1,j+1)*dx(i-1))/dxdy(i,j)
                       case ('23')  ! x ノイマン y 周期  u(i,j-1)=u(i,j)-f(i,j-1)*dy(j-1)
                          ! u(1,1)=u(nx-1,1)=u(nx-1,2)-f(nx-1,1)*dy(1)
                          tmp_b=0.25*bt(i,j)*(psi(i+1,j+1) +psi(nx-1,2)-bnd(nx-1,1)*dy(1) -psi(i+1,j-1)-psi(nx-1,j+1))/dxdy(i,j)
                       case ('31')  ! x 周期 y 固定  u(i,j-1)=u(i,ny-1)
                          tmp_b=0.25*bt(i,j)*(psi(i+1,j+1)+psi(i-1,j-1) -psi(i+1,ny-1)-psi(i-1,j+1))/dxdy(i,j)
                       case ('32')  ! x 周期 y ノイマン  u(i-1,j)=u(i,j)-f(i-1,j)*dx(i-1)
                          ! u(1,1)=u(1,ny-1)=u(2,ny-1)-f(1,ny-1)*dx(1)
                          ! u(i,j-1)=u(i,ny-1)
                          tmp_b=0.25*bt(i,j)*(psi(i+1,j+1) +psi(i,ny-1)-bnd(i-1,ny-1)*dx(i-1) -psi(i+1,ny-1)-psi(i,j+1)+bnd(i-1,j+1)*dx(i-1))/dxdy(i,j)
                       case ('33')  ! x, y 周期  u(1,1)=u(nx-1,ny-1)
                          tmp_b=0.25*bt(i,j)*(psi(i+1,j+1)+psi(nx-1,ny-1) -psi(i+1,ny-1)-psi(nx-1,j+1))/dxdy(i,j)

                       end select
                    else
                       tmp_b=0.0
                    end if

                    tmp=(tmp+tmp_adm+tmp_cem+tmp_b)/ac(i,j)

                 else  ! 左上隅  adp, cem は境界に関係なく求められる.

                    tmp=adp(i,j)*psi(i+1,j)+cem(i,j)*psi(i,j-1) -rho(i,j)

                    select case (bound(3:3))  ! x 軸上端条件 -> cep 項の計算
                    case ('1')
                       tmp_cep=cep(i,j)*psi(i,j+1)

                    case ('2')  ! u(i,j+1)=u(i,j)+f(i,j+1)*dy(j+1)
                       tmp_cep=cep(i,j)*bnd(i,j+1)*dy(j+1)

                    case ('3')  ! u(i,ny)=u(i,2)
                       tmp_cep=cep(i,j)*psi(i,2)

                    end select

                    select case (bound(2:2))
                    case ('1')
                       tmp_adm=adm(i,j)*psi(i-1,j)

                    case ('2')  ! u(i-1,j)=u(i,j)-f(i-1,j)*dx(i-1)
                       tmp_adm=-adm(i,j)*bnd(i-1,j)*dx(i-1)

                    case ('3')  ! u(1,j)=u(nx-1,j)
                       tmp_adm=adm(i,j)*psi(nx-1,j)

                    end select

                    if(signb==1)then
                       select case (bound(2:3))
                       case ('11')  ! 両方固定境界
                          tmp_b=0.25*bt(i,j)*(psi(i+1,j+1)+psi(i-1,j-1) -psi(i+1,j-1)-psi(i-1,j+1))/dxdy(i,j)
                       case ('12')  ! x ノイマン y 固定
                          ! u(i,j+1)=u(i,j)+f(i,j+1)*dy(j+1), i-1 は固定強制.
                          tmp_b=0.25*bt(i,j)*(psi(i+1,j)+bnd(i+1,j+1)*dy(j+1) +psi(i-1,j-1) -psi(i+1,j-1)-psi(i-1,j+1))/dxdy(i,j)
                       case ('13')  ! x 周期 y 固定  u(i,ny)=u(i,2), i-1 は固定強制.
                          tmp_b=0.25*bt(i,j)*(psi(i+1,2)+psi(i-1,j-1) -psi(i+1,j-1)-psi(i-1,j+1))/dxdy(i,j)
                       case ('21')  ! x 固定 y ノイマン  u(i-1,j)=u(i,j)-f(i-1,j)*dx(i-1)
                          ! j+1 は固定強制.
                          tmp_b=0.25*bt(i,j)*(psi(i+1,j+1) +psi(i,j-1)-bnd(i-1,j-1)*dx(i-1) -psi(i+1,j-1)-psi(i-1,j+1))/dxdy(i,j)
                       case ('22')  ! x, y ノイマン  u(i-1,j)=u(i,j)-f(i-1,j)*dx(i-1)
                          ! u(i-1,j+1)=u(i,j)+0.5*(-f(i-1,j+1)*dx(i-1)+f(i-1,j+1)*dy(j+1))
                          ! u(i,j+1)=u(i,j)+f(i,j+1)*dy(j+1)
                          tmp_b=0.25*bt(i,j)*(psi(i+1,j)+bnd(i+1,j+1)*dy(j+1) +psi(i,j-1)-bnd(i-1,j-1)*dx(i-1) -psi(i+1,j-1) +0.5*(bnd(i-1,j+1)*dx(i-1)-bnd(i-1,j+1)*dy(j+1)))/dxdy(i,j)
                       case ('23')  ! x 周期 y ノイマン  u(i,ny)=u(i,2)
                          ! u(1,ny)=u(1,2)=u(2,2)-f(1,2)*dx(1)
                          ! u(i-1,j)=u(i,j)-f(i-1,j)*dx(i-1)
                          tmp_b=0.25*bt(i,j)*(psi(i+1,2) +psi(i,j-1)-bnd(i-1,j-1)*dx(i-1) -psi(i+1,j-1)-psi(2,2)+bnd(1,2)*dx(1))/dxdy(i,j)
                       case ('31')  ! x 固定 y 周期  u(1,j)=u(nx-1,j),  j+1 は固定強制
                          tmp_b=0.25*bt(i,j)*(psi(i+1,j+1)+psi(nx-1,j-1) -psi(i+1,j-1)-psi(i-1,j+1))/dxdy(i,j)
                       case ('32')  ! x ノイマン y 周期  u(i,j+1)=u(i,j)+f(i,j+1)*dy(j+1)
                          ! u(1,ny)=u(nx-1,ny)=u(nx-1,ny-1)+f(nx-1,ny)*dy(ny)
                          ! u(1,j)=u(nx-1,j)
                          tmp_b=0.25*bt(i,j)*(psi(i+1,j)+bnd(i+1,j+1)*dy(j+1) +psi(nx-1,j-1)-psi(i+1,j-1) -psi(nx-1,ny-1)-bnd(nx-1,ny)*dy(ny))/dxdy(i,j)
                       case ('33')  ! x, y 周期  u(1,ny)=u(nx-1,2)
                          tmp_b=0.25*bt(i,j)*(psi(i+1,2)+psi(nx-1,j-1) -psi(i+1,j-1)-psi(nx-1,2))/dxdy(i,j)

                       end select
                    else
                       tmp_b=0.0
                    end if

                    tmp=(tmp+tmp_adm+tmp_cep+tmp_b)/ac(i,j)


                 end if

              end if

           end if

           if(i==nx-1)then  ! 右境界での扱い. ここでは, adp 以外の項は共通計算可能.
              if(j/=2.or.j/=ny-1)then  ! 隅境界以外
                 tmp=cep(i,j)*psi(i,j+1) +adm(i,j)*psi(i-1,j)+cem(i,j)*psi(i,j-1)-rho(i,j)

                 select case (bound(4:4))  ! y 軸右端
                 case ('1')
                    tmp_adp=adp(i,j)*psi(i+1,j)

                    if(signb==1)then
                       tmp_b=0.25*bt(i,j)*(psi(i+1,j+1)+psi(i-1,j-1) -psi(i+1,j-1)-psi(i-1,j+1))/dxdy(i,j)
                    else
                       tmp_b=0.0
                    end if

                 case ('2')  ! u(i+1,j)=u(i,j)+f(i+1,j)*dx(i+1)
                    tmp_adp=adp(i,j)*bnd(i+1,j)*dx(i+1)

                    if(signb==1)then
                       tmp_b=0.25*bt(i,j)*(psi(i,j+1)+bnd(i+1,j+1)*dx(i+1) +psi(i-1,j-1)-psi(i,j-1)-bnd(i+1,j-1)*dx(i+1) -psi(i-1,j+1))/dxdy(i,j)
                    else
                       tmp_b=0.0
                    end if

                 case ('3')  ! u(nx,j)=u(2,j)
                    tmp_adp=adp(i,j)*psi(2,j)

                    if(signb==1)then
                       tmp_b=0.25*bt(i,j)*(psi(2,j+1)+psi(i-1,j-1) -psi(2,j-1)-psi(i-1,j+1))/dxdy(i,j)
                    else
                       tmp_b=0.0
                    end if

                 end select

                 tmp=(tmp+tmp_adp+tmp_b)/ac(i,j)

              else

                 if(j==2)then  ! 右下隅  adm, cep は境界に関係なく求められる.

                    tmp=adm(i,j)*psi(i-1,j)+cep(i,j)*psi(i,j+1) -rho(i,j)

                    select case (bound(1:1))  ! x 軸下端条件 -> cem 項の計算
                    case ('1')
                       tmp_cem=cem(i,j)*psi(i,j-1)

                    case ('2')  ! u(i,j-1)=u(i,j)-f(i,j-1)*dy(j-1)
                       tmp_cem=-cem(i,j)*bnd(i,j-1)*dy(j-1)

                    case ('3')  ! u(i,1)=u(i,ny-1)
                       tmp_cem=cem(i,j)*psi(i,ny-1)

                    end select

                    select case (bound(4:4))
                    case ('1')
                       tmp_adp=adp(i,j)*psi(i-1,j)

                    case ('2')  ! u(i+1,j)=u(i,j)+f(i+1,j)*dx(i+1)
                       tmp_adp=adp(i,j)*bnd(i+1,j)*dx(i+1)

                    case ('3')  ! u(nx,j)=u(2,j)
                       tmp_adp=adp(i,j)*psi(2,j)

                    end select

                    if(signb==1)then
                       select case (bound(1:1)//bound(4:4))
                       case ('11')  ! 両方固定境界
                          tmp_b=0.25*bt(i,j)*(psi(i+1,j+1)+psi(i-1,j-1) -psi(i+1,j-1)-psi(i-1,j+1))/dxdy(i,j)
                       case ('12')  ! x 固定 y ノイマン
                          ! u(i+1,j)=u(i,j)+f(i+1,j)*dx(i+1), j-1 は固定強制.
                          tmp_b=0.25*bt(i,j)*(psi(i,j+1)+bnd(i+1,j+1)*dx(i+1) +psi(i-1,j-1)-psi(i+1,j-1)-psi(i,j+1))/dxdy(i,j)
                       case ('13')  ! x 固定 y 周期  u(nx,j)=u(2,j)
                          tmp_b=0.25*bt(i,j)*(psi(2,j+1)+psi(i-1,j-1) -psi(i+1,j-1)-psi(i-1,j+1))/dxdy(i,j)
                       case ('21')  ! x ノイマン y 固定  u(i,j-1)=u(i,j)-f(i,j-1)*dy(j-1)
                          tmp_b=0.25*bt(i,j)*(psi(i+1,j+1)+psi(i-1,j-1) -psi(i+1,j)+bnd(i+1,j-1)*dy(j-1)-psi(i-1,j+1))/dxdy(i,j)
                       case ('22')  ! x, y ノイマン  u(i+1,j)=u(i,j)+f(i+1,j)*dx(i+1)
                          ! u(i+1,j-1)=u(i,j)+0.5*(f(i+1,j)*dx(i+1)-f(i,j-1)*dy(j-1))
                          ! u(i,j-1)=u(i,j)-f(i,j-1)*dy(j-1)
                          tmp_b=0.25*bt(i,j)*(psi(i,j+1)+bnd(i+1,j+1)*dx(i+1) -psi(i-1,j)+bnd(i-1,j-1)*dy(j-1) -0.5*(bnd(i+1,j)*dx(i+1)-bnd(i,j-1)*dy(j-1)) -psi(i-1,j+1))/dxdy(i,j)
                       case ('23')  ! x ノイマン y 周期  u(i,j-1)=u(i,j)-f(i,j-1)*dy(j-1)
                          ! u(nx,1)=u(2,1)=u(2,2)-f(2,1)*dy(1)
                          tmp_b=0.25*bt(i,j)*(psi(2,j+1) +psi(i-1,j)-bnd(i-1,j-1)*dy(j-1) -psi(2,2)+bnd(2,1)*dy(1)-psi(i-1,j+1))/dxdy(i,j)
                       case ('31')  ! x 周期 y 固定  u(i,j-1)=u(i,ny-1)
                          tmp_b=0.25*bt(i,j)*(psi(i+1,j+1)+psi(i-1,ny-1) -psi(i+1,j-1)-psi(i-1,j+1))/dxdy(i,j)
                       case ('32')  ! x 周期 y ノイマン  u(i-1,j)=u(i,j)-f(i-1,j)*dx(i-1)
                          ! u(1,1)=u(1,ny-1)=u(2,ny-1)-f(1,ny-1)*dx(1)
                          ! u(i,j-1)=u(i,ny-1)
                          tmp_b=0.25*bt(i,j)*(psi(i+1,j+1) +psi(i,ny-1)-bnd(i-1,ny-1)*dx(i-1) -psi(i+1,ny-1)-psi(i,j+1)+bnd(i-1,j+1)*dx(i-1))/dxdy(i,j)
                       case ('33')  ! x, y 周期  u(nx,1)=u(2,ny-1)
                          tmp_b=0.25*bt(i,j)*(psi(2,j+1)+psi(i-1,ny-1) -psi(2,ny-1)-psi(i-1,j+1))/dxdy(i,j)

                      end select
                    else
                       tmp_b=0.0
                    end if

                    tmp=(tmp+tmp_adp+tmp_cem+tmp_b)/ac(i,j)

                 else  ! 右上隅  adm, cem は境界に関係なく求められる.

                    tmp=adm(i,j)*psi(i-1,j)+cem(i,j)*psi(i,j-1) -rho(i,j)

                    select case (bound(3:3))  ! x 軸上端条件 -> cep 項の計算
                    case ('1')
                       tmp_cep=cep(i,j)*psi(i,j+1)

                    case ('2')  ! u(i,j+1)=u(i,j)+f(i,j+1)*dy(j+1)
                       tmp_cep=cep(i,j)*bnd(i,j+1)*dy(j+1)

                    case ('3')  ! u(i,ny)=u(i,2)
                       tmp_cep=cep(i,j)*psi(i,2)

                    end select

                    select case (bound(4:4))
                    case ('1')
                       tmp_adp=adp(i,j)*psi(i-1,j)

                    case ('2')  ! u(i+1,j)=u(i,j)+f(i+1,j)*dx(i+1)
                       tmp_adp=adp(i,j)*bnd(i+1,j)*dx(i+1)

                    case ('3')  ! u(nx,j)=u(2,j)
                       tmp_adp=adp(i,j)*psi(2,j)

                    end select

                    if(signb==1)then
                       select case (bound(3:4))
                       case ('11')  ! 両方固定境界
                          tmp_b=0.25*bt(i,j)*(psi(i+1,j+1)+psi(i-1,j-1) -psi(i+1,j-1)-psi(i-1,j+1))/dxdy(i,j)
                       case ('12')  ! x 固定 y ノイマン  u(i+1,j)=u(i,j)+f(i+1,j)*dx(i+1)
                          ! j+1 は固定強制.
                          tmp_b=0.25*bt(i,j)*(psi(i+1,j+1)+psi(i-1,j-1) -psi(i+1,j)-bnd(i+1,j-1)*dx(i+1)-psi(i-1,j+1))/dxdy(i,j)
                       case ('13')  ! x 固定 y 周期  u(nx,j)=u(2,j),  j+1 は固定強制
                          tmp_b=0.25*bt(i,j)*(psi(i+1,j+1)+psi(i-1,j-1) -psi(2,j-1)-psi(i-1,j+1))/dxdy(i,j)
                       case ('21')  ! x ノイマン y 固定
                          ! u(i,j+1)=u(i,j)+f(i,j+1)*dy(j+1), i+1 は固定強制.
                          tmp_b=0.25*bt(i,j)*(psi(i+1,j+1)+psi(i-1,j-1) -psi(i+1,j-1)-psi(i-1,j)-bnd(i-1,j+1)*dy(j+1))/dxdy(i,j)
                       case ('22')  ! x, y ノイマン  u(i+1,j)=u(i,j)+f(i+1,j)*dx(i+1)
                          ! u(i+1,j+1)=u(i,j)+0.5*(f(i+1,j)*dx(i+1)+f(i,j+1)*dy(j+1))
                          ! u(i,j+1)=u(i,j)+f(i,j+1)*dy(j+1)
                          tmp_b=0.25*bt(i,j)*(0.5*(bnd(i+1,j)*dx(i+1) +bnd(i,j+1)*dy(j+1))+psi(i-1,j-1) -psi(i-1,j)-bnd(i-1,j+1)*dy(j+1) -psi(i,j-1)-bnd(i+1,j-1)*dx(i+1))/dxdy(i,j)
                       case ('23')  ! x ノイマン y 周期  u(i,j+1)=u(i,j)+f(i,j+1)*dy(j+1)
                          ! u(nx,ny)=u(2,ny)=u(2,ny-1)+f(2,ny)*dy(ny)
                          ! u(nx,j)=u(2,j)
                          tmp_b=0.25*bt(i,j)*(psi(2,ny-1)+bnd(2,ny)*dy(ny) +psi(i-1,j-1)-psi(2,j-1) -psi(i-1,j)-bnd(i-1,j+1)*dy(j+1))/dxdy(i,j)
                       case ('31')  ! x 周期 y 固定  u(i,ny)=u(i,2), i+1 は固定強制.
                          tmp_b=0.25*bt(i,j)*(psi(i+1,j+1)+psi(i-1,j-1) -psi(i+1,j-1)-psi(i-1,2))/dxdy(i,j)
                       case ('32')  ! x 周期 y ノイマン  u(i,ny)=u(i,2)
                          ! u(nx,ny)=u(nx,2)=u(nx-1,2)+f(nx,2)*dx(nx)
                          ! u(i+1,j)=u(i,j)+f(i+1,j)*dx(i+1)
                          tmp_b=0.25*bt(i,j)*(psi(nx-1,2)+bnd(nx,2)*dx(nx) +psi(i-1,j-1)-psi(i,j-1)-bnd(i+1,j-1)*dx(i+1) -psi(i-1,2))/dxdy(i,j)
                       case ('33')  ! x, y 周期  u(nx,ny)=u(2,2)
                          tmp_b=0.25*bt(i,j)*(psi(2,2)+psi(i-1,j-1) -psi(2,j-1)-psi(i-1,2))/dxdy(i,j)

                       end select
                    else
                       tmp_b=0.0
                    end if

                    tmp=(tmp+tmp_adp+tmp_cep+tmp_b)/ac(i,j)


                 end if

              end if
           end if

           if(j==2)then  ! 下境界での扱い. ここでは, cem 以外の項は共通計算できる.
              if(i/=2.or.i/=nx-1)then  ! 隅境界以外
                 tmp=adp(i,j)*psi(i+1,j)+cep(i,j)*psi(i,j+1) +adm(i,j)*psi(i-1,j)-rho(i,j)

                 select case (bound(1:1))  ! x 軸下端
                 case ('1')
                    tmp_cem=cem(i,j)*psi(i,j-1)

                    if(signb==1)then
                       tmp_b=0.25*bt(i,j)*(psi(i+1,j+1)+psi(i-1,j-1) -psi(i+1,j-1)-psi(i-1,j+1))/dxdy(i,j)
                    else
                       tmp_b=0.0
                    end if

                 case ('2')  ! u(i,j-1)=u(i,j)-f(i,j-1)*dy(j-1)
                    tmp_cem=-cem(i,j)*bnd(i,j-1)*dy(j-1)
            
                    if(signb==1)then
                       tmp_b=0.25*bt(i,j)*(psi(i+1,j+1) +psi(i-1,j)-bnd(i-1,j-1)*dy(j-1) -psi(i+1,j)+bnd(i+1,j-1)*dy(j-1)-psi(i-1,j+1))/dxdy(i,j)
                    else
                       tmp_b=0.0
                    end if

                 case ('3')  ! u(i,1)=u(i,ny-1)
                    tmp_cem=cem(i,j)*psi(i,ny-1)

                    if(signb==1)then
                       tmp_b=0.25*bt(i,j)*(psi(i+1,j+1)+psi(i-1,ny-1) -psi(i+1,ny-1)-psi(i-1,j+1))/dxdy(i,j)
                    else
                       tmp_b=0.0
                    end if

                 end select

                 tmp=(tmp+tmp_cem+tmp_b)/ac(i,j)

              end if
           end if

           if(j==ny-1)then  ! 上境界での扱い. ここでは, cep 以外の項は共通計算可能.
              if(i/=2.or.i/=nx-1)then  ! 隅境界以外
                 tmp=adp(i,j)*psi(i+1,j) +adm(i,j)*psi(i-1,j)+cem(i,j)*psi(i,j-1)-rho(i,j)

                 select case (bound(3:3))  ! x 軸上端
                 case ('1')
                    tmp_cep=cep(i,j)*psi(i,j+1)

                    if(signb==1)then
                       tmp_b=0.25*bt(i,j)*(psi(i+1,j+1)+psi(i-1,j-1) -psi(i+1,j-1)-psi(i-1,j+1))/dxdy(i,j)
                    else
                       tmp_b=0.0
                    end if

                 case ('2')  ! u(i,j+1)=u(i,j)+f(i,j+1)*dy(j+1)
                    tmp_cep=cep(i,j)*bnd(i,j+1)*dy(j+1)

                    if(signb==1)then
                       tmp_b=0.25*bt(i,j)*(psi(i+1,j)+bnd(i+1,j+1)*dy(j+1) +psi(i-1,j-i)-psi(i+1,j-1) -psi(i-1,j)-bnd(i-1,j+1)*dy(j+1))/dxdy(i,j)
                    else
                       tmp_b=0.0
                    end if

                 case ('3')  ! u(i,ny)=u(i,2)
                    tmp_cep=cep(i,j)*psi(i,2)

                    if(signb==1)then
                       tmp_b=0.25*bt(i,j)*(psi(i+1,2)+psi(i-1,j-1) -psi(i+1,j-1)-psi(i-1,2))/dxdy(i,j)
                    else
                       tmp_b=0.0
                    end if

                 end select

                 tmp=(tmp+tmp_cep+tmp_b)/ac(i,j)

              end if

           end if
        end if
!!-- ここまでは if 文のインデントを調節しない.
!        end if
!
!        tmp=tmp/ac(i,j)
!-- 誤差の計算 ---
!        if(psi(i,j)==0.0)then
!           if(tmp==0.0)then
!              err=0.0
!           else
!              err=abs(tmp-psi(i,j))/abs(tmp)
!           end if
!        else
!           if(abs(psi(i,j))<1.0.and.abs(tmp)<1.0)then
!              err=abs(tmp-psi(i,j))
!           else
!              err=abs(tmp-psi(i,j))/abs(psi(i,j))
!           end if
!        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
!-- 境界の設定
!-- x 下端境界
  select case (bound(1:1))  ! bound = 1 は境界での更新がないので, 場合分けしない.
  case ('2')
     do i=2,nx-1
        psi(i,1)=psi(i,2)-bnd(i,1)*dy(1)
     end do

  case ('3')
     do i=2,nx-1
        psi(i,1)=psi(i,ny-1)
     end do
  end select

!-- y 左端境界
  select case (bound(2:2))  ! bound = 1 は境界での更新がないので, 場合分けしない.
  case ('2')
     do j=2,ny-1
        psi(1,j)=psi(2,j)-bnd(1,j)*dx(1)
     end do

  case ('3')
     do j=2,ny-1
        psi(1,j)=psi(nx-1,j)
     end do
  end select

!-- x 上端境界
  select case (bound(3:3))  ! bound = 1 は境界での更新がないので, 場合分けしない.
  case ('2')
     do i=2,nx-1
        psi(i,ny)=psi(i,ny-1)+bnd(i,ny)*dy(ny)
     end do

  case ('3')
     do i=2,nx-1
        psi(i,ny)=psi(i,2)
     end do
  end select

!-- y 右端境界
  select case (bound(4:4))  ! bound = 1 は境界での更新がないので, 場合分けしない.
  case ('2')
     do j=2,ny-1
        psi(nx,j)=psi(nx-1,j)+bnd(nx,j)*dx(nx)
     end do

  case ('3')
     do j=2,ny-1
        psi(nx,j)=psi(2,j)
     end do
  end select

!-- 4 端  ! 以下, 固定境界がどちらかにあるなら, その値に強制 ('11', '12,', '13', '21', '31')
!-- また, bound(1:2) の値が分かれば, bound(3:4) の値も一意に決まる.
!-- 片端のみ周期境界なら, 周期境界で値を更新し, そのあとにフラックスを流す順番にする.
  select case (bound(1:2))
  case ('23')  ! この場合は y 軸両端で周期境界なので, (1,1), (nx,1) の値がここで決まる.
     psi(1,1)=psi(nx-1,2)-bnd(nx-1,1)*dy(1)
     psi(nx,1)=psi(2,2)-bnd(2,1)*dy(1)

     if(bound(3:3)=='2')then  ! x 軸上端が 1 ならすでに強制されているので, 何もしない.
        psi(1,ny)=psi(nx-1,ny-1)+bnd(nx-1,ny)*dy(ny)
        psi(nx,ny)=psi(2,ny-1)+bnd(2,ny)*dy(ny)
     end if

  case ('32')  ! この場合は x 軸両端で周期境界なので, (1,1), (1,ny) の値がここで決まる.
     psi(1,1)=psi(2,ny-1)-bnd(1,ny-1)*dx(1)
     psi(1,ny)=psi(2,2)-bnd(1,2)*dx(1)

     if(bound(4:4)=='2')then  ! y 軸右端が 1 ならすでに強制されているので, 何もしない.
        psi(nx,1)=psi(nx-1,ny-1)+bnd(nx,ny-1)*dx(nx)
        psi(nx,ny)=psi(nx-1,2)+bnd(nx,2)*dx(nx)
     end if

  case ('33')  ! この場合は全領域周期境界
     psi(1,1)=psi(nx-1,ny-1)
     psi(nx,1)=psi(2,ny-1)
     psi(1,ny)=psi(nx-1,2)
     psi(nx,ny)=psi(2,2)

  case default  ! あとは, どの隅においても, '22' でなければ, 固定強制される.
     if(bound(1:2)=='22')then
        psi(1,1)=psi(2,2)-0.5*(bnd(1,2)*dy(1)+bnd(2,1)*dx(1))
     end if

     if(bound(2:3)=='22')then
        psi(1,ny)=psi(2,ny-1)+0.5*(bnd(1,ny-1)*dy(ny)-bnd(2,ny)*dx(1))
     end if

     if(bound(3:4)=='22')then
        psi(nx,ny)=psi(nx-1,ny-1)+0.5*(bnd(nx,ny-1)*dy(ny)+bnd(nx-1,ny)*dx(nx))
     end if

     if(bound(1:1)=='2'.and.bound(4:4)=='2')then
        psi(nx,1)=psi(nx-1,2)+0.5*(-bnd(nx,2)*dy(1)+bnd(nx-1,1)*dx(nx))
     end if

  end select
  end do
  

end subroutine Poisson_GauSei
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
: 各微分項の係数

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

[Source]

subroutine Poisson_Jacobi(x, y, rho, eps, boundary, psi, bound_opt, a, b, c, d, e)
  ! ヤコビ法によるポアソン方程式の求積(開発中openmp)
  ! 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(inout) :: psi(size(x),size(y))  ! ポアソン方程式の解
  integer :: i, j, k, l, m, n
  integer :: nx  ! x 方向の配列要素
  integer :: ny  ! y 方向の配列要素
  integer :: signb, signd, signe  ! 各係数を計算するかどうか
  real :: err, err_max
  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
  real, dimension(size(x),size(y)) :: adp, adm, cep, cem, ac
  character(4) :: bound
  real, dimension(size(x),size(y)) :: tmp, tmp_b, tmp_adm, tmp_adp, tmp_cem, tmp_cep

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

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

  psi = 0.0

!-- 周期境界の設定確認.
!-- 周期境界なので, 両端とも 3 が設定されていないといけない.
  if(bound(1:1)=='3')then
     if(bound(3:3)/='3')then
        write(*,*) "### ERROR ###"
        write(*,*) "if bound = 3, bound(1:1)==bound(3:3). STOP."
        stop
     end if
  end if

  if(bound(3:3)=='3')then
     if(bound(1:1)/='3')then
        write(*,*) "### ERROR ###"
        write(*,*) "if bound = 3, bound(1:1)==bound(3:3). STOP."
        stop
     end if
  end if

  if(bound(2:2)=='3')then
     if(bound(4:4)/='3')then
        write(*,*) "### ERROR ###"
        write(*,*) "if bound = 3, bound(2:2)==bound(4:4). STOP."
        stop
     end if
  end if

  if(bound(4:4)=='3')then
     if(bound(2:2)/='3')then
        write(*,*) "### ERROR ###"
        write(*,*) "if bound = 3, bound(2:2)==bound(4:4). STOP."
        stop
     end if
  end if

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

  if(present(c))then
     do j=1,ny
        do i=1,nx
           ct(i,j)=c(i,j)
        end do
     end do
  else
     do j=1,ny
        do i=1,nx
           ct(i,j)=1.0
        end do
     end do
  end if

  if(present(b))then
     do j=1,ny
        do i=1,nx
           bt(i,j)=b(i,j)
        end do
     end do
     signb=1
  else
     signb=0
  end if

  if(present(d))then
     do j=1,ny
        do i=1,nx
           dt(i,j)=d(i,j)
        end do
     end do
     signd=1
  else
     signd=0
  end if

  if(present(e))then
     do j=1,ny
        do i=1,nx
           et(i,j)=e(i,j)
        end do
     end do
     signe=1
  else
     signe=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

!-- ポアソン係数の計算

  if(signd==0)then  ! 付加項 d がついていないとき
     do j=2,ny-1
        do i=2,nx-1
           adp(i,j)=at(i,j)/(dx2(i))
           adm(i,j)=at(i,j)/(dx2(i))
        end do
     end do
  else
     do j=2,ny-1
        do i=2,nx-1
           adp(i,j)=at(i,j)/(dx2(i))+0.5*dt(i,j)/dx(i)
           adm(i,j)=at(i,j)/(dx2(i))-0.5*dt(i,j)/dx(i)
        end do
     end do
  end if

  if(signe==0)then  ! 付加項 e がついていないとき
     do j=2,ny-1
        do i=2,nx-1
           cep(i,j)=ct(i,j)/(dy2(j))
           cem(i,j)=ct(i,j)/(dy2(j))
        end do
     end do
  else
     do j=2,ny-1
        do i=2,nx-1
           cep(i,j)=ct(i,j)/(dy2(j))+0.5*et(i,j)/dy(j)
           cem(i,j)=ct(i,j)/(dy2(j))-0.5*et(i,j)/dy(j)
        end do
     end do
  end if

!-- 最高次数係数 ac の計算 (境界条件によって評価式が変わる.)
!-- 実際にポアソン係数が使用されるのは, 境界を除いた 1 回り内側なので,
!-- 計算量削減のため, ループをそのようにしておく.
  do j=2,ny-1  ! 係数 ac については, 境界の内側で境界条件によって式が異なる.
     do i=2,nx-1  ! ここでは, 固定端条件での値を代入しており, それ以外の場合は以下の処理で上書きする.
        ac(i,j)=2.0*(at(i,j)/dx2(i)+ct(i,j)/dy2(j))
     end do
  end do

!-- 以下, 境界の内側での係数 ac の計算
!-- 本計算では, 固定境界があるなら, 隅領域はすべてその値で固定するようにしており,
!-- 固定境界値は反復計算の前に与えられる値を用いるので, 境界条件 1 では評価式は変わらない.
!-- また, 周期境界条件では, 両端が周期境界でなければならないから, bound の組み合わせは
!-- 制限される.
!-- 実際に係数が変わるのは, 条件 2 のときだけであるので, 2 について評価式を計算する.
!-- (3 は係数ではなく, ポアソンの式だけが変化する.)
!-- 以下では, 境界の 4 隅の内側も計算している.
!-- これで 2 重勘定されるのは, 隅の両側ともが 2 の条件のときだけなので,
!-- この方法を用いた方が場合分けが少なくて済む.
!-- 2 重勘定したときは, その半分を引く.

  if(bound(1:1)=='2')then  ! x 下端がノイマン条件 : psi(i,1)=psi(i,2)-f(i,1)*dy(1)
     do i=2,nx-1
        ac(i,2)=ac(i,2)-cem(i,2)
     end do
  end if

  if(bound(2:2)=='2')then  ! y 左端がノイマン条件 : psi(1,j)=psi(2,j)-f(1,j)*dx(1)
     do j=2,ny-1
        ac(2,j)=ac(2,j)-adm(2,j)
     end do
  end if

  if(bound(3:3)=='2')then  ! x 上端がノイマン条件 : psi(i,ny)=psi(i,ny-1)+f(i,ny)*dy(ny)
     do i=2,nx-1
        ac(i,ny-1)=ac(i,ny-1)-cep(i,ny-1)
     end do
  end if

  if(bound(4:4)=='2')then  ! y 右端がノイマン条件 : psi(nx,j)=psi(nx-1,j)+f(nx,j)*dx(nx)
     do j=2,ny-1
        ac(nx-1,j)=ac(nx-1,j)-adp(nx-1,j)
     end do
  end if

!-- 内側 4 隅での 2 重勘定の解消  bt についてはここでのみ計算されるので, 足し合わせておく.

  if(bound(1:2)=='22')then  ! 左下隅の解消
     ac(2,2)=ac(2,2)+(cem(2,2)+adm(2,2))*0.5-0.25*bt(2,2)/dxdy(2,2)
  end if

  if(bound(2:3)=='22')then  ! 左上隅の解消
     ac(2,ny-1)=ac(2,ny-1)+(cep(2,ny-1)+adm(2,ny-1))*0.5+0.25*bt(2,ny-1)/dxdy(2,ny-1)
  end if

  if(bound(3:4)=='22')then  ! 右上隅の解消
     ac(nx-1,ny-1)=ac(nx-1,ny-1)+(cep(nx-1,ny-1)+adp(nx-1,ny-1))*0.5 -0.25*bt(nx-1,ny-1)/dxdy(nx-1,ny-1)
  end if

  if(bound(1:1)=='2'.and.bound(4:4)=='2')then  ! 右下隅の解消
     ac(nx-1,2)=ac(nx-1,2)+(cem(nx-1,2)+adp(nx-1,2))*0.5+0.25*bt(nx-1,2)/dxdy(nx-1,2)
  end if

!-- 境界値の設定
  if(present(bound_opt))then
     do j=1,ny
        do i=1,nx
           bnd(i,j)=bound_opt(i,j)
        end do
     end do
  else
     do j=1,ny
        do i=1,nx
           bnd(i,j)=0.0
        end do
     end do
  end if

!-- 境界条件の確認 (+ 固定境界なら, その値を代入する.)
  if(bound(1:1)=='1')then
     do i=1,nx
        psi(i,1)=bnd(i,1)
     end do
  end if

  if(bound(2:2)=='1')then
     do j=1,ny
        psi(1,j)=bnd(1,j)
     end do
  end if

  if(bound(3:3)=='1')then
     do i=1,nx
        psi(i,ny)=bnd(i,ny)
     end do
  end if

  if(bound(4:4)=='1')then
     do j=1,ny
        psi(1,j)=bnd(1,j)
     end do
  end if

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

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

        if(i/=2.and.i/=nx-1.and.j/=2.and.j/=nx-1)then
           ! 境界の 1 つ内側以外
           tmp(i,j)=adp(i,j)*psi(i+1,j)+cep(i,j)*psi(i,j+1) +adm(i,j)*psi(i-1,j)+cem(i,j)*psi(i,j-1) -rho(i,j)

           if(signb==1)then
              tmp(i,j)=tmp(i,j)+0.25*bt(i,j)*(psi(i+1,j+1)+psi(i-1,j-1) -psi(i+1,j-1)-psi(i-1,j+1))/dxdy(i,j)
           end if

           tmp(i,j)=tmp(i,j)/ac(i,j)

        else

!-- 以下で境界の内側計算
!-- コードを見やすくするため, adm, adp, cem, cep の項ごとに case 分けする.

           if(i==2)then  ! 左境界での扱い. ここでは, adm 以外の項は共通計算できる.
              if(j/=2.or.j/=ny-1)then  ! 隅境界以外
                 tmp(i,j)=adp(i,j)*psi(i+1,j)+cep(i,j)*psi(i,j+1) +cem(i,j)*psi(i,j-1)-rho(i,j)

                 select case (bound(2:2))  ! y 軸左端
                 case ('1')
                    tmp_adm(i,j)=adm(i,j)*psi(i-1,j)

                    if(signb==1)then
                       tmp_b(i,j)=0.25*bt(i,j)*(psi(i+1,j+1)+psi(i-1,j-1) -psi(i+1,j-1)-psi(i-1,j+1))/dxdy(i,j)
                    else
                       tmp_b(i,j)=0.0
                    end if

                 case ('2')  ! u(i-1,j)=u(i,j)-f(i-1,j)*dx(i-1)
                    tmp_adm(i,j)=-adm(i,j)*bnd(i-1,j)*dx(i-1)
            
                    if(signb==1)then
                       tmp_b(i,j)=0.25*bt(i,j)*(psi(i+1,j+1) +psi(i,j-1)-bnd(i-1,j-1)*dx(i-1) -psi(i+1,j-1)-psi(i,j+1)+bnd(i-1,j+1)*dx(i-1))/dxdy(i,j)
                    else
                       tmp_b(i,j)=0.0
                    end if

                 case ('3')  ! u(1,j)=u(nx-1,j)
                    tmp_adm(i,j)=adm(i,j)*psi(nx-1,j)

                    if(signb==1)then
                       tmp_b(i,j)=0.25*bt(i,j)*(psi(i+1,j+1)+psi(nx-1,j-1) -psi(i+1,j-1)-psi(nx-1,j+1))/dxdy(i,j)
                    else
                       tmp_b(i,j)=0.0
                    end if

                 end select

                 tmp(i,j)=(tmp(i,j)+tmp_adm(i,j)+tmp_b(i,j))/ac(i,j)

              else

                 if(j==2)then  ! 左下隅  ! adp, cep は境界条件に依存しない.
                    tmp(i,j)=adp(i,j)*psi(i+1,j)+cep(i,j)*psi(i,j+1) -rho(i,j)

                    select case (bound(1:1))  ! x 軸下端条件 -> cem 項の計算
                    case ('1')
                       tmp_cem(i,j)=cem(i,j)*psi(i,j-1)

                    case ('2')  ! u(i,j-1)=u(i,j)-f(i,j-1)*dy(j-1)
                       tmp_cem(i,j)=-cem(i,j)*bnd(i,j-1)*dy(j-1)

                    case ('3')  ! u(i,1)=u(i,ny-1)
                       tmp_cem(i,j)=cem(i,j)*psi(i,ny-1)

                    end select

                    select case (bound(2:2))
                    case ('1')
                       tmp_adm(i,j)=adm(i,j)*psi(i-1,j)

                    case ('2')  ! u(i-1,j)=u(i,j)-f(i-1,j)*dx(i-1)
                       tmp_adm(i,j)=-adm(i,j)*bnd(i-1,j)*dx(i-1)

                    case ('3')  ! u(1,j)=u(nx-1,j)
                       tmp_adm(i,j)=adm(i,j)*psi(nx-1,j)

                    end select

                    if(signb==1)then
                       select case (bound(1:2))
                       case ('11')  ! 両方固定境界
                          tmp_b(i,j)=0.25*bt(i,j)*(psi(i+1,j+1)+psi(i-1,j-1) -psi(i+1,j-1)-psi(i-1,j+1))/dxdy(i,j)
                       case ('12')  ! x 固定 y ノイマン
                          ! u(i-1,j)=u(i,j)-f(i-1,j)*dx(i-1), j-1 は固定強制.
                          tmp_b(i,j)=0.25*bt(i,j)*(psi(i+1,j+1)+psi(i-1,j-1) -psi(i+1,j-1)-psi(i,j+1)+bnd(i-1,j+1)*dx(i-1))/dxdy(i,j)
                       case ('13')  ! x 固定 y 周期  u(1,j)=u(nx-1,j)
                          tmp_b(i,j)=0.25*bt(i,j)*(psi(i+1,j+1)+psi(i-1,j-1) -psi(i+1,j-1)-psi(nx-1,j+1))/dxdy(i,j)
                       case ('21')  ! x ノイマン y 固定  u(i,j-1)=u(i,j)-f(i,j-1)*dy(j-1)
                          tmp_b(i,j)=0.25*bt(i,j)*(psi(i+1,j+1)+psi(i-1,j-1) -psi(i+1,j)+bnd(i+1,j-1)*dy(j-1)-psi(i-1,j+1))/dxdy(i,j)
                       case ('22')  ! x, y ノイマン  u(i-1,j)=u(i,j)-f(i-1,j)*dx(i-1)
                          ! u(i-1,j-1)=u(i,j)-0.5*(f(i-1,j)*dx(i-1)+f(i,j-1)*dy(j-1))
                          ! u(i,j-1)=u(i,j)-f(i,j-1)*dy(j-1)
                          tmp_b(i,j)=0.25*bt(i,j)*(psi(i+1,j+1) -0.5*(bnd(i-1,j)*dx(i-1)+bnd(i,j-1)*dy(j-1)) -psi(i+1,j)+bnd(i+1,j-1)*dy(j-1) -psi(i,j+1)+bnd(i-1,j+1)*dx(i-1))/dxdy(i,j)
                       case ('23')  ! x ノイマン y 周期  u(i,j-1)=u(i,j)-f(i,j-1)*dy(j-1)
                          ! u(1,1)=u(nx-1,1)=u(nx-1,2)-f(nx-1,1)*dy(1)
                          tmp_b(i,j)=0.25*bt(i,j)*(psi(i+1,j+1) +psi(nx-1,2)-bnd(nx-1,1)*dy(1) -psi(i+1,j-1)-psi(nx-1,j+1))/dxdy(i,j)
                       case ('31')  ! x 周期 y 固定  u(i,j-1)=u(i,ny-1)
                          tmp_b(i,j)=0.25*bt(i,j)*(psi(i+1,j+1)+psi(i-1,j-1) -psi(i+1,ny-1)-psi(i-1,j+1))/dxdy(i,j)
                       case ('32')  ! x 周期 y ノイマン  u(i-1,j)=u(i,j)-f(i-1,j)*dx(i-1)
                          ! u(1,1)=u(1,ny-1)=u(2,ny-1)-f(1,ny-1)*dx(1)
                          ! u(i,j-1)=u(i,ny-1)
                          tmp_b(i,j)=0.25*bt(i,j)*(psi(i+1,j+1) +psi(i,ny-1)-bnd(i-1,ny-1)*dx(i-1) -psi(i+1,ny-1)-psi(i,j+1)+bnd(i-1,j+1)*dx(i-1))/dxdy(i,j)
                       case ('33')  ! x, y 周期  u(1,1)=u(nx-1,ny-1)
                          tmp_b(i,j)=0.25*bt(i,j)*(psi(i+1,j+1)+psi(nx-1,ny-1) -psi(i+1,ny-1)-psi(nx-1,j+1))/dxdy(i,j)

                       end select
                    else
                       tmp_b(i,j)=0.0
                    end if

                    tmp(i,j)=(tmp(i,j)+tmp_adm(i,j)+tmp_cem(i,j)+tmp_b(i,j))/ac(i,j)

                 else  ! 左上隅  adp, cem は境界に関係なく求められる.

                    tmp(i,j)=adp(i,j)*psi(i+1,j)+cem(i,j)*psi(i,j-1) -rho(i,j)

                    select case (bound(3:3))  ! x 軸上端条件 -> cep 項の計算
                    case ('1')
                       tmp_cep(i,j)=cep(i,j)*psi(i,j+1)

                    case ('2')  ! u(i,j+1)=u(i,j)+f(i,j+1)*dy(j+1)
                       tmp_cep(i,j)=cep(i,j)*bnd(i,j+1)*dy(j+1)

                    case ('3')  ! u(i,ny)=u(i,2)
                       tmp_cep(i,j)=cep(i,j)*psi(i,2)

                    end select

                    select case (bound(2:2))
                    case ('1')
                       tmp_adm(i,j)=adm(i,j)*psi(i-1,j)

                    case ('2')  ! u(i-1,j)=u(i,j)-f(i-1,j)*dx(i-1)
                       tmp_adm(i,j)=-adm(i,j)*bnd(i-1,j)*dx(i-1)

                    case ('3')  ! u(1,j)=u(nx-1,j)
                       tmp_adm(i,j)=adm(i,j)*psi(nx-1,j)

                    end select

                    if(signb==1)then
                       select case (bound(2:3))
                       case ('11')  ! 両方固定境界
                          tmp_b(i,j)=0.25*bt(i,j)*(psi(i+1,j+1)+psi(i-1,j-1) -psi(i+1,j-1)-psi(i-1,j+1))/dxdy(i,j)
                       case ('12')  ! x ノイマン y 固定
                          ! u(i,j+1)=u(i,j)+f(i,j+1)*dy(j+1), i-1 は固定強制.
                          tmp_b(i,j)=0.25*bt(i,j)*(psi(i+1,j)+bnd(i+1,j+1)*dy(j+1) +psi(i-1,j-1) -psi(i+1,j-1)-psi(i-1,j+1))/dxdy(i,j)
                       case ('13')  ! x 周期 y 固定  u(i,ny)=u(i,2), i-1 は固定強制.
                          tmp_b(i,j)=0.25*bt(i,j)*(psi(i+1,2)+psi(i-1,j-1) -psi(i+1,j-1)-psi(i-1,j+1))/dxdy(i,j)
                       case ('21')  ! x 固定 y ノイマン  u(i-1,j)=u(i,j)-f(i-1,j)*dx(i-1)
                          ! j+1 は固定強制.
                          tmp_b(i,j)=0.25*bt(i,j)*(psi(i+1,j+1) +psi(i,j-1)-bnd(i-1,j-1)*dx(i-1) -psi(i+1,j-1)-psi(i-1,j+1))/dxdy(i,j)
                       case ('22')  ! x, y ノイマン  u(i-1,j)=u(i,j)-f(i-1,j)*dx(i-1)
                          ! u(i-1,j+1)=u(i,j)+0.5*(-f(i-1,j+1)*dx(i-1)+f(i-1,j+1)*dy(j+1))
                          ! u(i,j+1)=u(i,j)+f(i,j+1)*dy(j+1)
                          tmp_b(i,j)=0.25*bt(i,j)*(psi(i+1,j)+bnd(i+1,j+1)*dy(j+1) +psi(i,j-1)-bnd(i-1,j-1)*dx(i-1) -psi(i+1,j-1) +0.5*(bnd(i-1,j+1)*dx(i-1)-bnd(i-1,j+1)*dy(j+1)))/dxdy(i,j)
                       case ('23')  ! x 周期 y ノイマン  u(i,ny)=u(i,2)
                          ! u(1,ny)=u(1,2)=u(2,2)-f(1,2)*dx(1)
                          ! u(i-1,j)=u(i,j)-f(i-1,j)*dx(i-1)
                          tmp_b(i,j)=0.25*bt(i,j)*(psi(i+1,2) +psi(i,j-1)-bnd(i-1,j-1)*dx(i-1) -psi(i+1,j-1)-psi(2,2)+bnd(1,2)*dx(1))/dxdy(i,j)
                       case ('31')  ! x 固定 y 周期  u(1,j)=u(nx-1,j),  j+1 は固定強制
                          tmp_b(i,j)=0.25*bt(i,j)*(psi(i+1,j+1)+psi(nx-1,j-1) -psi(i+1,j-1)-psi(i-1,j+1))/dxdy(i,j)
                       case ('32')  ! x ノイマン y 周期  u(i,j+1)=u(i,j)+f(i,j+1)*dy(j+1)
                          ! u(1,ny)=u(nx-1,ny)=u(nx-1,ny-1)+f(nx-1,ny)*dy(ny)
                          ! u(1,j)=u(nx-1,j)
                          tmp_b(i,j)=0.25*bt(i,j)*(psi(i+1,j)+bnd(i+1,j+1)*dy(j+1) +psi(nx-1,j-1)-psi(i+1,j-1) -psi(nx-1,ny-1)-bnd(nx-1,ny)*dy(ny))/dxdy(i,j)
                       case ('33')  ! x, y 周期  u(1,ny)=u(nx-1,2)
                          tmp_b(i,j)=0.25*bt(i,j)*(psi(i+1,2)+psi(nx-1,j-1) -psi(i+1,j-1)-psi(nx-1,2))/dxdy(i,j)

                       end select
                    else
                       tmp_b(i,j)=0.0
                    end if

                    tmp(i,j)=(tmp(i,j)+tmp_adm(i,j)+tmp_cep(i,j)+tmp_b(i,j))/ac(i,j)


                 end if

              end if

           end if

           if(i==nx-1)then  ! 右境界での扱い. ここでは, adp 以外の項は共通計算可能.
              if(j/=2.or.j/=ny-1)then  ! 隅境界以外
                 tmp(i,j)=cep(i,j)*psi(i,j+1) +adm(i,j)*psi(i-1,j)+cem(i,j)*psi(i,j-1)-rho(i,j)

                 select case (bound(4:4))  ! y 軸右端
                 case ('1')
                    tmp_adp(i,j)=adp(i,j)*psi(i+1,j)

                    if(signb==1)then
                       tmp_b(i,j)=0.25*bt(i,j)*(psi(i+1,j+1)+psi(i-1,j-1) -psi(i+1,j-1)-psi(i-1,j+1))/dxdy(i,j)
                    else
                       tmp_b(i,j)=0.0
                    end if

                 case ('2')  ! u(i+1,j)=u(i,j)+f(i+1,j)*dx(i+1)
                    tmp_adp(i,j)=adp(i,j)*bnd(i+1,j)*dx(i+1)

                    if(signb==1)then
                       tmp_b(i,j)=0.25*bt(i,j)*(psi(i,j+1)+bnd(i+1,j+1)*dx(i+1) +psi(i-1,j-1)-psi(i,j-1)-bnd(i+1,j-1)*dx(i+1) -psi(i-1,j+1))/dxdy(i,j)
                    else
                       tmp_b(i,j)=0.0
                    end if

                 case ('3')  ! u(nx,j)=u(2,j)
                    tmp_adp(i,j)=adp(i,j)*psi(2,j)

                    if(signb==1)then
                       tmp_b(i,j)=0.25*bt(i,j)*(psi(2,j+1)+psi(i-1,j-1) -psi(2,j-1)-psi(i-1,j+1))/dxdy(i,j)
                    else
                       tmp_b(i,j)=0.0
                    end if

                 end select

                 tmp(i,j)=(tmp(i,j)+tmp_adp(i,j)+tmp_b(i,j))/ac(i,j)

              else

                 if(j==2)then  ! 右下隅  adm, cep は境界に関係なく求められる.

                    tmp(i,j)=adm(i,j)*psi(i-1,j)+cep(i,j)*psi(i,j+1) -rho(i,j)

                    select case (bound(1:1))  ! x 軸下端条件 -> cem 項の計算
                    case ('1')
                       tmp_cem(i,j)=cem(i,j)*psi(i,j-1)

                    case ('2')  ! u(i,j-1)=u(i,j)-f(i,j-1)*dy(j-1)
                       tmp_cem(i,j)=-cem(i,j)*bnd(i,j-1)*dy(j-1)

                    case ('3')  ! u(i,1)=u(i,ny-1)
                       tmp_cem(i,j)=cem(i,j)*psi(i,ny-1)

                    end select

                    select case (bound(4:4))
                    case ('1')
                       tmp_adp(i,j)=adp(i,j)*psi(i-1,j)

                    case ('2')  ! u(i+1,j)=u(i,j)+f(i+1,j)*dx(i+1)
                       tmp_adp(i,j)=adp(i,j)*bnd(i+1,j)*dx(i+1)

                    case ('3')  ! u(nx,j)=u(2,j)
                       tmp_adp(i,j)=adp(i,j)*psi(2,j)

                    end select

                    if(signb==1)then
                       select case (bound(1:1)//bound(4:4))
                       case ('11')  ! 両方固定境界
                          tmp_b(i,j)=0.25*bt(i,j)*(psi(i+1,j+1)+psi(i-1,j-1) -psi(i+1,j-1)-psi(i-1,j+1))/dxdy(i,j)
                       case ('12')  ! x 固定 y ノイマン
                          ! u(i+1,j)=u(i,j)+f(i+1,j)*dx(i+1), j-1 は固定強制.
                          tmp_b(i,j)=0.25*bt(i,j)*(psi(i,j+1)+bnd(i+1,j+1)*dx(i+1) +psi(i-1,j-1)-psi(i+1,j-1)-psi(i,j+1))/dxdy(i,j)
                       case ('13')  ! x 固定 y 周期  u(nx,j)=u(2,j)
                          tmp_b(i,j)=0.25*bt(i,j)*(psi(2,j+1)+psi(i-1,j-1) -psi(i+1,j-1)-psi(i-1,j+1))/dxdy(i,j)
                       case ('21')  ! x ノイマン y 固定  u(i,j-1)=u(i,j)-f(i,j-1)*dy(j-1)
                          tmp_b(i,j)=0.25*bt(i,j)*(psi(i+1,j+1)+psi(i-1,j-1) -psi(i+1,j)+bnd(i+1,j-1)*dy(j-1)-psi(i-1,j+1))/dxdy(i,j)
                       case ('22')  ! x, y ノイマン  u(i+1,j)=u(i,j)+f(i+1,j)*dx(i+1)
                          ! u(i+1,j-1)=u(i,j)+0.5*(f(i+1,j)*dx(i+1)-f(i,j-1)*dy(j-1))
                          ! u(i,j-1)=u(i,j)-f(i,j-1)*dy(j-1)
                          tmp_b(i,j)=0.25*bt(i,j)*(psi(i,j+1)+bnd(i+1,j+1)*dx(i+1) -psi(i-1,j)+bnd(i-1,j-1)*dy(j-1) -0.5*(bnd(i+1,j)*dx(i+1)-bnd(i,j-1)*dy(j-1)) -psi(i-1,j+1))/dxdy(i,j)
                       case ('23')  ! x ノイマン y 周期  u(i,j-1)=u(i,j)-f(i,j-1)*dy(j-1)
                          ! u(nx,1)=u(2,1)=u(2,2)-f(2,1)*dy(1)
                          tmp_b(i,j)=0.25*bt(i,j)*(psi(2,j+1) +psi(i-1,j)-bnd(i-1,j-1)*dy(j-1) -psi(2,2)+bnd(2,1)*dy(1)-psi(i-1,j+1))/dxdy(i,j)
                       case ('31')  ! x 周期 y 固定  u(i,j-1)=u(i,ny-1)
                          tmp_b(i,j)=0.25*bt(i,j)*(psi(i+1,j+1)+psi(i-1,ny-1) -psi(i+1,j-1)-psi(i-1,j+1))/dxdy(i,j)
                       case ('32')  ! x 周期 y ノイマン  u(i-1,j)=u(i,j)-f(i-1,j)*dx(i-1)
                          ! u(1,1)=u(1,ny-1)=u(2,ny-1)-f(1,ny-1)*dx(1)
                          ! u(i,j-1)=u(i,ny-1)
                          tmp_b(i,j)=0.25*bt(i,j)*(psi(i+1,j+1) +psi(i,ny-1)-bnd(i-1,ny-1)*dx(i-1) -psi(i+1,ny-1)-psi(i,j+1)+bnd(i-1,j+1)*dx(i-1))/dxdy(i,j)
                       case ('33')  ! x, y 周期  u(nx,1)=u(2,ny-1)
                          tmp_b(i,j)=0.25*bt(i,j)*(psi(2,j+1)+psi(i-1,ny-1) -psi(2,ny-1)-psi(i-1,j+1))/dxdy(i,j)

                      end select
                    else
                       tmp_b(i,j)=0.0
                    end if

                    tmp(i,j)=(tmp(i,j)+tmp_adp(i,j)+tmp_cem(i,j)+tmp_b(i,j))/ac(i,j)

                 else  ! 右上隅  adm, cem は境界に関係なく求められる.

                    tmp(i,j)=adm(i,j)*psi(i-1,j)+cem(i,j)*psi(i,j-1) -rho(i,j)

                    select case (bound(3:3))  ! x 軸上端条件 -> cep 項の計算
                    case ('1')
                       tmp_cep(i,j)=cep(i,j)*psi(i,j+1)

                    case ('2')  ! u(i,j+1)=u(i,j)+f(i,j+1)*dy(j+1)
                       tmp_cep(i,j)=cep(i,j)*bnd(i,j+1)*dy(j+1)

                    case ('3')  ! u(i,ny)=u(i,2)
                       tmp_cep(i,j)=cep(i,j)*psi(i,2)

                    end select

                    select case (bound(4:4))
                    case ('1')
                       tmp_adp(i,j)=adp(i,j)*psi(i-1,j)

                    case ('2')  ! u(i+1,j)=u(i,j)+f(i+1,j)*dx(i+1)
                       tmp_adp(i,j)=adp(i,j)*bnd(i+1,j)*dx(i+1)

                    case ('3')  ! u(nx,j)=u(2,j)
                       tmp_adp(i,j)=adp(i,j)*psi(2,j)

                    end select

                    if(signb==1)then
                       select case (bound(3:4))
                       case ('11')  ! 両方固定境界
                          tmp_b(i,j)=0.25*bt(i,j)*(psi(i+1,j+1)+psi(i-1,j-1) -psi(i+1,j-1)-psi(i-1,j+1))/dxdy(i,j)
                       case ('12')  ! x 固定 y ノイマン  u(i+1,j)=u(i,j)+f(i+1,j)*dx(i+1)
                          ! j+1 は固定強制.
                          tmp_b(i,j)=0.25*bt(i,j)*(psi(i+1,j+1)+psi(i-1,j-1) -psi(i+1,j)-bnd(i+1,j-1)*dx(i+1)-psi(i-1,j+1))/dxdy(i,j)
                       case ('13')  ! x 固定 y 周期  u(nx,j)=u(2,j),  j+1 は固定強制
                          tmp_b(i,j)=0.25*bt(i,j)*(psi(i+1,j+1)+psi(i-1,j-1) -psi(2,j-1)-psi(i-1,j+1))/dxdy(i,j)
                       case ('21')  ! x ノイマン y 固定
                          ! u(i,j+1)=u(i,j)+f(i,j+1)*dy(j+1), i+1 は固定強制.
                          tmp_b(i,j)=0.25*bt(i,j)*(psi(i+1,j+1)+psi(i-1,j-1) -psi(i+1,j-1)-psi(i-1,j)-bnd(i-1,j+1)*dy(j+1))/dxdy(i,j)
                       case ('22')  ! x, y ノイマン  u(i+1,j)=u(i,j)+f(i+1,j)*dx(i+1)
                          ! u(i+1,j+1)=u(i,j)+0.5*(f(i+1,j)*dx(i+1)+f(i,j+1)*dy(j+1))
                          ! u(i,j+1)=u(i,j)+f(i,j+1)*dy(j+1)
                          tmp_b(i,j)=0.25*bt(i,j)*(0.5*(bnd(i+1,j)*dx(i+1) +bnd(i,j+1)*dy(j+1))+psi(i-1,j-1) -psi(i-1,j)-bnd(i-1,j+1)*dy(j+1) -psi(i,j-1)-bnd(i+1,j-1)*dx(i+1))/dxdy(i,j)
                       case ('23')  ! x ノイマン y 周期  u(i,j+1)=u(i,j)+f(i,j+1)*dy(j+1)
                          ! u(nx,ny)=u(2,ny)=u(2,ny-1)+f(2,ny)*dy(ny)
                          ! u(nx,j)=u(2,j)
                          tmp_b(i,j)=0.25*bt(i,j)*(psi(2,ny-1)+bnd(2,ny)*dy(ny) +psi(i-1,j-1)-psi(2,j-1) -psi(i-1,j)-bnd(i-1,j+1)*dy(j+1))/dxdy(i,j)
                       case ('31')  ! x 周期 y 固定  u(i,ny)=u(i,2), i+1 は固定強制.
                          tmp_b(i,j)=0.25*bt(i,j)*(psi(i+1,j+1)+psi(i-1,j-1) -psi(i+1,j-1)-psi(i-1,2))/dxdy(i,j)
                       case ('32')  ! x 周期 y ノイマン  u(i,ny)=u(i,2)
                          ! u(nx,ny)=u(nx,2)=u(nx-1,2)+f(nx,2)*dx(nx)
                          ! u(i+1,j)=u(i,j)+f(i+1,j)*dx(i+1)
                          tmp_b(i,j)=0.25*bt(i,j)*(psi(nx-1,2)+bnd(nx,2)*dx(nx) +psi(i-1,j-1)-psi(i,j-1)-bnd(i+1,j-1)*dx(i+1) -psi(i-1,2))/dxdy(i,j)
                       case ('33')  ! x, y 周期  u(nx,ny)=u(2,2)
                          tmp_b(i,j)=0.25*bt(i,j)*(psi(2,2)+psi(i-1,j-1) -psi(2,j-1)-psi(i-1,2))/dxdy(i,j)

                       end select
                    else
                       tmp_b(i,j)=0.0
                    end if

                    tmp(i,j)=(tmp(i,j)+tmp_adp(i,j)+tmp_cep(i,j)+tmp_b(i,j))/ac(i,j)


                 end if

              end if
           end if

           if(j==2)then  ! 下境界での扱い. ここでは, cem 以外の項は共通計算できる.
              if(i/=2.or.i/=nx-1)then  ! 隅境界以外
                 tmp(i,j)=adp(i,j)*psi(i+1,j)+cep(i,j)*psi(i,j+1) +adm(i,j)*psi(i-1,j)-rho(i,j)

                 select case (bound(1:1))  ! x 軸下端
                 case ('1')
                    tmp_cem(i,j)=cem(i,j)*psi(i,j-1)

                    if(signb==1)then
                       tmp_b(i,j)=0.25*bt(i,j)*(psi(i+1,j+1)+psi(i-1,j-1) -psi(i+1,j-1)-psi(i-1,j+1))/dxdy(i,j)
                    else
                       tmp_b(i,j)=0.0
                    end if

                 case ('2')  ! u(i,j-1)=u(i,j)-f(i,j-1)*dy(j-1)
                    tmp_cem(i,j)=-cem(i,j)*bnd(i,j-1)*dy(j-1)
            
                    if(signb==1)then
                       tmp_b(i,j)=0.25*bt(i,j)*(psi(i+1,j+1) +psi(i-1,j)-bnd(i-1,j-1)*dy(j-1) -psi(i+1,j)+bnd(i+1,j-1)*dy(j-1)-psi(i-1,j+1))/dxdy(i,j)
                    else
                       tmp_b(i,j)=0.0
                    end if

                 case ('3')  ! u(i,1)=u(i,ny-1)
                    tmp_cem(i,j)=cem(i,j)*psi(i,ny-1)

                    if(signb==1)then
                       tmp_b(i,j)=0.25*bt(i,j)*(psi(i+1,j+1)+psi(i-1,ny-1) -psi(i+1,ny-1)-psi(i-1,j+1))/dxdy(i,j)
                    else
                       tmp_b(i,j)=0.0
                    end if

                 end select

                 tmp(i,j)=(tmp(i,j)+tmp_cem(i,j)+tmp_b(i,j))/ac(i,j)

              end if
           end if

           if(j==ny-1)then  ! 上境界での扱い. ここでは, cep 以外の項は共通計算可能.
              if(i/=2.or.i/=nx-1)then  ! 隅境界以外
                 tmp(i,j)=adp(i,j)*psi(i+1,j) +adm(i,j)*psi(i-1,j)+cem(i,j)*psi(i,j-1)-rho(i,j)

                 select case (bound(3:3))  ! x 軸上端
                 case ('1')
                    tmp_cep(i,j)=cep(i,j)*psi(i,j+1)

                    if(signb==1)then
                       tmp_b(i,j)=0.25*bt(i,j)*(psi(i+1,j+1)+psi(i-1,j-1) -psi(i+1,j-1)-psi(i-1,j+1))/dxdy(i,j)
                    else
                       tmp_b(i,j)=0.0
                    end if

                 case ('2')  ! u(i,j+1)=u(i,j)+f(i,j+1)*dy(j+1)
                    tmp_cep(i,j)=cep(i,j)*bnd(i,j+1)*dy(j+1)

                    if(signb==1)then
                       tmp_b(i,j)=0.25*bt(i,j)*(psi(i+1,j)+bnd(i+1,j+1)*dy(j+1) +psi(i-1,j-i)-psi(i+1,j-1) -psi(i-1,j)-bnd(i-1,j+1)*dy(j+1))/dxdy(i,j)
                    else
                       tmp_b(i,j)=0.0
                    end if

                 case ('3')  ! u(i,ny)=u(i,2)
                    tmp_cep(i,j)=cep(i,j)*psi(i,2)

                    if(signb==1)then
                       tmp_b(i,j)=0.25*bt(i,j)*(psi(i+1,2)+psi(i-1,j-1) -psi(i+1,j-1)-psi(i-1,2))/dxdy(i,j)
                    else
                       tmp_b(i,j)=0.0
                    end if

                 end select

                 tmp(i,j)=(tmp(i,j)+tmp_cep(i,j)+tmp_b(i,j))/ac(i,j)

              end if

           end if
        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
!        if(psi(i,j)==0.0)then
!           err=abs(tmp(i,j)-psi(i,j))/abs(tmp(i,j))
!        else
!           err=abs(tmp(i,j)-psi(i,j))/abs(psi(i,j))
!        end if
        err=abs(tmp(i,j)-psi(i,j))

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

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

!-- 境界の設定
!-- x 下端境界
  select case (bound(1:1))  ! bound = 1 は境界での更新がないので, 場合分けしない.
  case ('2')
     do i=2,nx-1
        psi(i,1)=psi(i,2)-bnd(i,1)*dy(1)
     end do

  case ('3')
     do i=2,nx-1
        psi(i,1)=psi(i,ny-1)
     end do
  end select

!-- y 左端境界
  select case (bound(2:2))  ! bound = 1 は境界での更新がないので, 場合分けしない.
  case ('2')
     do j=2,ny-1
        psi(1,j)=psi(2,j)-bnd(1,j)*dx(1)
     end do

  case ('3')
     do j=2,ny-1
        psi(1,j)=psi(nx-1,j)
     end do
  end select

!-- x 上端境界
  select case (bound(3:3))  ! bound = 1 は境界での更新がないので, 場合分けしない.
  case ('2')
     do i=2,nx-1
        psi(i,ny)=psi(i,ny-1)+bnd(i,ny)*dy(ny)
     end do

  case ('3')
     do i=2,nx-1
        psi(i,ny)=psi(i,2)
     end do
  end select

!-- y 右端境界
  select case (bound(4:4))  ! bound = 1 は境界での更新がないので, 場合分けしない.
  case ('2')
     do j=2,ny-1
        psi(nx,j)=psi(nx-1,j)+bnd(nx,j)*dx(nx)
     end do

  case ('3')
     do j=2,ny-1
        psi(nx,j)=psi(2,j)
     end do
  end select

!-- 4 端  ! 以下, 固定境界がどちらかにあるなら, その値に強制 ('11', '12,', '13', '21', '31')
!-- また, bound(1:2) の値が分かれば, bound(3:4) の値も一意に決まる.
!-- 片端のみ周期境界なら, 周期境界で値を更新し, そのあとにフラックスを流す順番にする.
  select case (bound(1:2))
  case ('23')  ! この場合は y 軸両端で周期境界なので, (1,1), (nx,1) の値がここで決まる.
     psi(1,1)=psi(nx-1,2)-bnd(nx-1,1)*dy(1)
     psi(nx,1)=psi(2,2)-bnd(2,1)*dy(1)

     if(bound(3:3)=='2')then  ! x 軸上端が 1 ならすでに強制されているので, 何もしない.
        psi(1,ny)=psi(nx-1,ny-1)+bnd(nx-1,ny)*dy(ny)
        psi(nx,ny)=psi(2,ny-1)+bnd(2,ny)*dy(ny)
     end if

  case ('32')  ! この場合は x 軸両端で周期境界なので, (1,1), (1,ny) の値がここで決まる.
     psi(1,1)=psi(2,ny-1)-bnd(1,ny-1)*dx(1)
     psi(1,ny)=psi(2,2)-bnd(1,2)*dx(1)

     if(bound(4:4)=='2')then  ! y 軸右端が 1 ならすでに強制されているので, 何もしない.
        psi(nx,1)=psi(nx-1,ny-1)+bnd(nx,ny-1)*dx(nx)
        psi(nx,ny)=psi(nx-1,2)+bnd(nx,2)*dx(nx)
     end if

  case ('33')  ! この場合は全領域周期境界
     psi(1,1)=psi(nx-1,ny-1)
     psi(nx,1)=psi(2,ny-1)
     psi(1,ny)=psi(nx-1,2)
     psi(nx,ny)=psi(2,2)

  case default  ! あとは, どの隅においても, '22' でなければ, 固定強制される.
     if(bound(1:2)=='22')then
        psi(1,1)=psi(2,2)-0.5*(bnd(1,2)*dy(1)+bnd(2,1)*dx(1))
     end if

     if(bound(2:3)=='22')then
        psi(1,ny)=psi(2,ny-1)+0.5*(bnd(1,ny-1)*dy(ny)-bnd(2,ny)*dx(1))
     end if

     if(bound(3:4)=='22')then
        psi(nx,ny)=psi(nx-1,ny-1)+0.5*(bnd(nx,ny-1)*dy(ny)+bnd(nx-1,ny)*dx(nx))
     end if

     if(bound(1:1)=='2'.and.bound(4:4)=='2')then
        psi(nx,1)=psi(nx-1,2)+0.5*(-bnd(nx,2)*dy(1)+bnd(nx-1,1)*dx(nx))
     end if

  end select
  end do
  
end subroutine Poisson_Jacobi