!---------------------------------------------------------------------
!     Copyright (C) GFD Dennou Club, 2006. All rights reserved.
!---------------------------------------------------------------------
!
!= 1 次元 (y 方向) 不等間隔交互格子 有限差分モデル用 境界条件モジュール
!
!* 履歴
!  * 2007/07/15 (小高正嗣) : 3D deepconv へ移植, dc_types を Use.
!  * 2006/06/02 (小高正嗣) : x_bc_module より作成
!

module y_bc_module
  != 1 次元 (y 方向) 不等間隔交互格子 有限差分モデル用 境界条件モジュール
  !
  !== 概要
  !
  ! y_bc_module は, 1 次元 (y 方向) 不等間隔交互格子を用いた有限差分法に
  ! 基づく数値モデルのための, 境界条件設定 Fortran 90 副プログラムを提供
  ! する.
  !
  ! このモジュールは y_module の下位モジュールである. 下請けモジュール
  ! として data_type, y_base_module モジュールを用いている. 
  !
  !
  !== 手続きの命名法
  !
  ! 境界条件を設定する副プログラムは
  !
  !   Boundary[境界条件を示す文字列]_(入力配列の次元情報)
  !
  ! のように命名されている. 境界条件を示す文字列は以下の 3 つである.
  !
  !   Sym  : 境界で対称
  !   Asym : 境界で反対称
  !   Cyc : 周期境界 
  ! 

  use dc_types,      only : DBKIND => DP
  use y_base_module, only : jm, ymargin, jmin, jmax
  implicit none

  private
  public :: BoundarySym_y, BoundaryAsym_y, BoundaryCyc_y
  public :: BoundarySym_q, BoundaryAsym_q, BoundaryCyc_q

  interface BoundaryCyc_y
    module procedure BoundaryCyc
  end interface 

  interface BoundaryCyc_q
    module procedure BoundaryCyc
  end interface 

  contains
!--------------------------------------------------------------------
    subroutine BoundarySym_y(y_Var)
      ! 1 次元 (y 方向) 配列に対称境界条件を適用する.

      real(DBKIND),intent(inout) :: y_Var(jmin:jmax) ! 入出力配列 
      integer                    :: jy               ! ループ添字

      ! 対称境界条件を適用する
      ! 
      !   ex.) jm=5, ymargin=2 の場合 (jmin=-1, jmax=7)
      !
      !   q    -1   0   1   2   3   4   5   6   7 
      !     |-:-|-:-|-:-|-:-|-:-|-:-|-:-|-:-|-:-|
      !   y  -1   0   1   2   3   4   5   6   7
      !     <------->                   <------->  
      !    y_Var(jmin:0)     と     y_Var(jm+1,jmax) の値を更新
      !     
      !   y_Var(0)  = y_Var(1)
      !   y_Var(-1) = y_Var(2)
      !   y_Var(6)  = y_Var(5)
      !   y_Var(7)  = y_Var(4)
      !
      do jy = 1, ymargin
        y_Var(1-jy)  = y_Var(jy)
        y_Var(jm+jy) = y_Var(jm+1-jy)
      end do

    end subroutine BoundarySym_y
!--------------------------------------------------------------------
    subroutine BoundaryAsym_y(y_Var)
      ! 1 次元 (y 方向) 配列に反対称境界条件を適用する.

      real(DBKIND),intent(inout) :: y_Var(jmin:jmax) ! 入出力配列 
      integer                    :: jy               ! ループ添字

      ! 反対称境界条件を適用する
      ! 
      !   ex.) jm=5, ymargin=2 の場合 (jmin=-1, jmax=7)
      !
      !   q    -1   0   1   2   3   4   5   6   7 
      !     |-:-|-:-|-:-|-:-|-:-|-:-|-:-|-:-|-:-|
      !   y  -1   0   1   2   3   4   5   6   7
      !     <------->                   <------->  
      !    y_Var(jmin:0)     と     y_Var(jm+1,jmax) の値を更新
      !     
      !   y_Var(0)  = - y_Var(1)
      !   y_Var(-1) = - y_Var(2)
      !   y_Var(6)  = - y_Var(5)
      !   y_Var(7)  = - y_Var(4)
      !
      do jy = 1, ymargin
        y_Var(1-jy)  = - y_Var(jy)
        y_Var(jm+jy) = - y_Var(jm+1-jy)
      end do

    end subroutine BoundaryAsym_y
!--------------------------------------------------------------------
    subroutine BoundarySym_q(q_Var)
      ! 1 次元 (y 方向) 配列に対称境界条件を適用する.

      real(DBKIND),intent(inout) :: q_Var(jmin:jmax) ! 入出力配列 
      integer                    :: jy               ! ループ添字

      ! 対称境界条件を適用する
      ! 
      !   ex.) jm=5, ymargin=2 の場合 (jmin=-1, jmax=7)
      !
      !   q    -1   0   1   2   3   4   5   6   7 
      !     |-:-|-:-|-:-|-:-|-:-|-:-|-:-|-:-|-:-|
      !   y  -1   0   1   2   3   4   5   6   7
      !     <------->                   <------->  
      !    q_Var(jmin:-1)    と     q_Var(jm+1:jmax) の値を更新
      !     
      !   q_Var(-1) = y_Var(1)
      !   q_Var(6)  = y_Var(4)
      !   q_Var(7)  = y_Var(3)
      !
      do jy = 1, ymargin-1
        q_Var(-jy)   = q_Var(jy)      
      end do

      do jy = 1, ymargin
        q_Var(jm+jy) = q_Var(jm+1-jy)
      end do

    end subroutine BoundarySym_q
!--------------------------------------------------------------------
    subroutine BoundaryAsym_q(q_Var)
      ! 1 次元 (y 方向) 配列に反対称境界条件を適用する.

      real(DBKIND),intent(inout) :: q_Var(jmin:jmax) ! 入出力配列 
      integer                    :: jy               ! ループ添字

      ! 反対称境界条件を適用する
      ! 
      !   ex.) jm=5, ymargin=2 の場合 (jmin=-1, jmax=7)
      !
      !   q    -1   0   1   2   3   4   5   6   7 
      !     |-:-|-:-|-:-|-:-|-:-|-:-|-:-|-:-|-:-|
      !   y  -1   0   1   2   3   4   5   6   7
      !     <------->                   <------->  
      !    q_Var(jmin:0)     と      q_Var(jm,jmax) の値を更新
      !     
      !   y_Var(0)  = 0
      !   y_Var(-1) = - y_Var(1)
      !   y_Var(5)  = 0
      !   y_Var(6)  = - y_Var(4)
      !   y_Var(7)  = - y_Var(3)
      !
      q_Var(0) = 0.0d0
      q_Var(jm) = 0.0d0
          
      do jy = 1, ymargin-1
        q_Var(-jy) = - q_Var(jy)
      end do

      do jy = 1, ymargin
        q_Var(jm+jy) = - q_Var(jm-jy)
      end do

    end subroutine BoundaryAsym_q
!--------------------------------------------------------------------
    subroutine BoundaryCyc(Var)
      ! 1 次元 (y 方向) 配列に周期境界条件を適用する.

      real(DBKIND),intent(inout) :: Var(jmin:jmax) ! 入出力配列 
      integer                    :: jy             ! ループ添字

      ! 周期界条件を適用する
      ! 
      !   ex.) jm=5, ymargin=2 の場合 (jmin=-1, jmax=7)
      !
      !   q    -1   0   1   2   3   4   5   6   7 
      !     |-:-|-:-|-:-|-:-|-:-|-:-|-:-|-:-|-:-|
      !   y  -1   0   1   2   3   4   5   6   7
      !     <------->                   <------->  
      !    y_Var(jmin:0)     と    y_Var(jm+1,jmax) の値を更新
      !    q_Var(jmin:0)     と    q_Var(jm+1,jmax) の値を更新
      !     
      !   y_Var(0)  = y_Var(5)
      !   y_Var(-1) = y_Var(4)
      !   y_Var(6)  = y_Var(1)
      !   y_Var(7)  = y_Var(2)
      !     
      !   q_Var(0)  = q_Var(5)
      !   q_Var(-1) = q_Var(4)
      !   q_Var(6)  = q_Var(1)
      !   q_Var(7)  = q_Var(2)
      !
      do jy = 1, ymargin
        Var(1-jy)  = Var(jm+1-jy)
        Var(jm+jy) = Var(jy)
      end do

    end subroutine BoundaryCyc
!--------------------------------------------------------------------
end module y_bc_module
