!---------------------------------------------------------------------
!     Copyright (C) GFD Dennou Club, 2004, 2006. All rights reserved.
!---------------------------------------------------------------------
!= Module Boundary
!
!   * Developer: SUGIYAMA Ko-ichiro 
!   * Version: $Id: boundary.f90,v 1.2 2005/04/22 15:01:52 sugiyama Exp $ 
!   * Tag Name: $Name:  $
!   * Change History: 
!
!== Overview 
!
!境界条件を与えるためのモジュール. 
!引数 type に指定された境界条件から, 「のり代」の部分の値を設定する
!
!== Error Handling
!
!== Known Bugs
!
!== Note
!
!水平鉛直に摩擦無し境界とする設定は, 移流のテスト計算用であることに注意.
!
!== Future Plans
!

module boundary

  !モジュール読み込み
  use gridset, only: Margin,        &! 境界のグリッド数
    &                DimXMin,       &! x 方向の配列の下限
    &                DimXMax,       &! x 方向の配列の上限
    &                DimZMin,       &! z 方向の配列の下限
    &                DimZMax,       &! z 方向の配列の上限
    &                RegXMin,       &! x 方向の物理領域の下限
    &                RegXMax,       &! x 方向の物理領域の上限
    &                RegZMin,       &! z 方向の物理領域の下限
    &                RegZMax         ! z 方向の物理領域の上限
     
  !暗黙の型宣言禁止
  implicit none

  !private 属性にする
  private
  
  !関数を public にする
  public xz_Boundary_Cyc_Sym
  public xz_Boundary_Cyc_AntiSym
  public pz_Boundary_Cyc_Sym
  public pz_Boundary_Cyc_AntiSym
  public xr_Boundary_Cyc_Sym
  public xr_Boundary_Cyc_AntiSym

contains

!!!---------------------------------------------------------------------!!!
  function Boundary_Cyc( Var )
    !
    ! x 方向に「周期境界条件」を適用する. 
    ! 格子点, 半格子点においても, 関数の形式は同じ
    !

    !暗黙の型宣言禁止
    implicit none
    
    !変数定義
    real(8)                 :: Boundary_Cyc(DimXMin:DimXMax)
    real(8), intent(inout)  :: Var(DimXMin:DimXMax)
    integer                 :: i


    Boundary_Cyc = Var

    Boundary_cyc(RegXMin) = Var(RegXMax)

    do i = 1, Margin
      Boundary_cyc(RegXMin - i) = Var(RegXMax - i)
      Boundary_cyc(RegXMax + i) = Var(RegXMin + i) 
    end do

  end function Boundary_Cyc


!!!---------------------------------------------------------------------!!!
  function z_Boundary_Sym( z_Var )
    !
    ! z 方向に半格子ずれた点に存在する変数に対し, 
    ! z 方向に「対称境界条件」を適用する. 
    !

    !暗黙の型宣言禁止
    implicit none
    
    !変数定義
    real(8)                 :: z_Boundary_Sym(DimZMin:DimZMax)
    real(8), intent(inout)  :: z_Var(DimZMin:DimZMax)
    integer                 :: k
  
    z_Boundary_Sym = z_Var

    do k = 0, Margin
      z_Boundary_Sym( RegZMin - k ) = z_Var( RegZMin + 1 + k )
    end do
    do k = 1, Margin
      z_Boundary_Sym( RegZMax + k ) = z_Var( RegZMax + 1 - k )
    end do
    
  end function z_Boundary_Sym


!!!---------------------------------------------------------------------!!!
  function z_Boundary_AntiSym( z_Var )
    !
    ! z 方向に半格子ずれた点に存在する変数に対し, 
    ! z 方向に「反対称境界条件」を適用する. 
    !

    !暗黙の型宣言禁止
    implicit none
    
    !変数定義
    real(8)                 :: z_Boundary_AntiSym(DimZMin:DimZMax)
    real(8), intent(inout)  :: z_Var(DimZMin:DimZMax)
    integer                 :: k
  
    z_Boundary_AntiSym = z_Var

    do k = 0, Margin
      z_Boundary_AntiSym( RegZMin - k ) = - z_Var( RegZMin + 1 + k )
    end do
    do k = 1, Margin
      z_Boundary_AntiSym( RegZMax + k ) = - z_Var( RegZMax + 1 - k )
    end do
    
  end function z_Boundary_AntiSym


!!!---------------------------------------------------------------------!!!
  function r_Boundary_Sym( r_Var )
    !
    ! z 方向の格子点上に存在する変数に対し, 
    ! z 方向に「対称境界条件」を適用する. 
    !
    
    !暗黙の型宣言禁止
    implicit none
    
    !変数定義
    real(8)                 :: r_Boundary_Sym(DimZMin:DimZMax)
    real(8), intent(inout)  :: r_Var(DimZMin:DimZMax)
    integer                 :: k
  
    r_Boundary_Sym = r_Var

    r_Boundary_Sym( RegZMin ) = 0.0d0
    r_Boundary_Sym( RegZMax ) = 0.0d0
    do k = 1, Margin
      r_Boundary_Sym( RegZMin - k ) = r_Var( RegZMin + k )
      r_Boundary_Sym( RegZMax + k ) = r_Var( RegZMax - k )
    end do
    
  end function r_Boundary_Sym


!!!---------------------------------------------------------------------!!!
  function r_Boundary_AntiSym( r_Var )
    !
    ! z 方向の格子点上に存在する変数に対し, 
    ! z 方向に「反対称境界条件」を適用する. 
    !
    
    !暗黙の型宣言禁止
    implicit none
    
    !変数定義
    real(8)                 :: r_Boundary_AntiSym(DimZMin:DimZMax)
    real(8), intent(inout)  :: r_Var(DimZMin:DimZMax)
    integer                 :: k
  
    r_Boundary_AntiSym = r_Var

    r_Boundary_AntiSym( RegZMin ) = 0.0d0
    r_Boundary_AntiSym( RegZMax ) = 0.0d0
    do k = 1, Margin
      r_Boundary_AntiSym( RegZMin - k ) = - r_Var( RegZMin + k )
      r_Boundary_AntiSym( RegZMax + k ) = - r_Var( RegZMax - k )
    end do
    
  end function r_Boundary_AntiSym


!!!---------------------------------------------------------------------!!!
  function xz_Boundary_Cyc_Sym(xz_Var)
    !
    ! x, z 方向に半格子ずれた点に存在する変数に対し, 
    ! x 方向に「周期境界条件」, z 方向に「対称境界条件」を適用する. 
    !

    !暗黙の型宣言禁止
    implicit none
    
    !変数定義
    real(8)        :: xz_Boundary_Cyc_Sym(DimXMin:DimXMax, DimZMin:DimZMax)
    real(8), intent(inout)      :: xz_var(DimXMin:DimXMax, DimZMin:DimZMax)
    integer                     :: i, k

    xz_Boundary_Cyc_Sym = xz_Var

    ! x 方向: 周期境界
    do k = RegZMin, RegZMax
      xz_Boundary_Cyc_Sym(:,k) = Boundary_Cyc( xz_Var(:,k) ) 
    end do

    !z 方向: 対称
    do i = DimXMin, DimXMax
      xz_Boundary_Cyc_Sym(i,:) = z_Boundary_Sym( xz_Boundary_Cyc_Sym(i,:) )
    end do

  end function xz_Boundary_Cyc_Sym


!!!---------------------------------------------------------------------!!!
  function xz_Boundary_Cyc_AntiSym(xz_Var)
    !
    ! x, z 方向に半格子ずれた点に存在する変数に対し, 
    ! x 方向に「周期境界条件」, z 方向に「反対称境界条件」を適用する. 
    !

    !暗黙の型宣言禁止
    implicit none
    
    !変数定義
    real(8)    :: xz_Boundary_Cyc_AntiSym(DimXMin:DimXMax, DimZMin:DimZMax)
    real(8), intent(inout)      :: xz_var(DimXMin:DimXMax, DimZMin:DimZMax)
    integer                     :: i, k
  
    xz_Boundary_Cyc_AntiSym = 0.0d0 
  
    ! x 方向: 周期境界 
    do k = RegZMin, RegZMax
      xz_Boundary_Cyc_AntiSym(:,k) = Boundary_Cyc( xz_Var(:,k) ) 
    end do

    !z 方向: 対称
    do i = DimXMin, DimXMax
      xz_Boundary_Cyc_AntiSym(i,:) =                         &
        & z_Boundary_AntiSym( xz_Boundary_Cyc_AntiSym(i,:))
    end do
    
  end function xz_Boundary_Cyc_AntiSym


!!!---------------------------------------------------------------------!!!
  function pz_Boundary_Cyc_Sym(pz_Var)
    !
    ! z 方向に半格子ずれた点に存在する変数に対し, 
    ! x 方向に「周期境界条件」, z 方向に「対称境界条件」を適用する. 
    !

    !暗黙の型宣言禁止
    implicit none
    
    !変数定義
    real(8)        :: pz_Boundary_Cyc_Sym(DimXMin:DimXMax, DimZMin:DimZMax)
    real(8), intent(inout)      :: pz_var(DimXMin:DimXMax, DimZMin:DimZMax)
    integer                     :: i, k

    pz_Boundary_Cyc_Sym = 0.0d0 
  
    ! x 方向: 周期境界
    do k = RegZMin, RegZMax
      pz_Boundary_Cyc_Sym(:,k) = Boundary_Cyc( pz_Var(:,k) ) 
    end do

    !z 方向: 対称
    do i = DimXMin, DimXMax
      pz_Boundary_Cyc_Sym(i,:) = z_Boundary_Sym( pz_Boundary_Cyc_Sym(i,:) ) 
    end do

  end function pz_Boundary_Cyc_Sym


!!!---------------------------------------------------------------------!!!
  function pz_Boundary_Cyc_AntiSym(pz_Var)
    !
    ! z 方向に半格子ずれた点に存在する変数に対し, 
    ! x 方向に「周期境界条件」, z 方向に「反対称境界条件」を適用する. 
    !

    !暗黙の型宣言禁止
    implicit none
    
    !変数定義
    real(8)    :: pz_Boundary_Cyc_AntiSym(DimXMin:DimXMax, DimZMin:DimZMax)
    real(8), intent(inout)      :: pz_var(DimXMin:DimXMax, DimZMin:DimZMax)
    integer                     :: i, k
  
    pz_Boundary_Cyc_AntiSym = 0.0d0 
  
    ! x 方向: 周期境界
    do k = RegZMin, RegZMax
      pz_Boundary_Cyc_AntiSym(:,k) = Boundary_Cyc( pz_Var(:,k) ) 
    end do

    !z 方向: 対称
    do i = DimXMin, DimXMax
      pz_Boundary_Cyc_AntiSym(i,:) =                         &
        & z_Boundary_AntiSym( pz_Boundary_Cyc_AntiSym(i,:) )
    end do
    
  end function pz_Boundary_Cyc_AntiSym


!!!---------------------------------------------------------------------!!!
  function xr_Boundary_Cyc_Sym( xr_Var )
    !
    ! x 方向に半格子ずれた点に存在する変数に対し, 
    ! x 方向に「周期境界条件」, z 方向に「対称境界条件」を適用する. 
    !

    !暗黙の型宣言禁止
    implicit none
    
    !変数定義
    real(8)        :: xr_Boundary_Cyc_Sym(DimXMin:DimXMax, DimZMin:DimZMax)
    real(8), intent(inout)      :: xr_var(DimXMin:DimXMax, DimZMin:DimZMax)
    integer                     :: i, k

    xr_Boundary_Cyc_Sym = 0.0d0 
  
    ! x 方向: 周期境界
    do k = RegZMin, RegZMax
      xr_Boundary_Cyc_Sym(:,k) = Boundary_Cyc( xr_Var(:,k) ) 
    end do

    !z 方向: 対称
    do i = DimXMin, DimXMax
      xr_Boundary_Cyc_Sym(i,:) = r_Boundary_Sym( xr_Boundary_Cyc_Sym(i,:) )
    end do
    
  end function xr_Boundary_Cyc_Sym


!!!---------------------------------------------------------------------!!!
  function xr_Boundary_Cyc_AntiSym( xr_Var )
    !
    ! x 方向に半格子ずれた点に存在する変数に対し, 
    ! x 方向に「周期境界条件」, z 方向に「反対称境界条件」を適用する. 
    !

    !暗黙の型宣言禁止
    implicit none
    
    !変数定義
    real(8)    :: xr_Boundary_Cyc_AntiSym(DimXMin:DimXMax, DimZMin:DimZMax)
    real(8), intent(inout)      :: xr_var(DimXMin:DimXMax, DimZMin:DimZMax)
    integer                     :: i, k
  
    xr_Boundary_Cyc_AntiSym = 0.0d0 
  
    ! x 方向: 周期境界
    do k = RegZMin, RegZMax
      xr_Boundary_Cyc_AntiSym(:,k) = Boundary_Cyc( xr_Var(:,k) ) 
    end do

    !z 方向: 反対称
    do i = DimXMin, DimXMax
      xr_Boundary_Cyc_AntiSym(i,:) =                         &
        & r_Boundary_AntiSym( xr_Boundary_Cyc_AntiSym(i,:) )
    end do

  end function xr_Boundary_Cyc_AntiSym

end module boundary
