!---------------------------------------------------------------------
!     Copyright (C) GFD Dennou Club, 2005. All rights reserved.
!---------------------------------------------------------------------
!= Module Damping
!
!   * Developer: SUGIYAMA Ko-ichiro
!   * Version: $Id: $
!   * Tag Name: $Name:  $
!   * Change History: 
!
!== Overview
!
!減衰率とその計算を行うためのパッケージ型モジュール
!  * 音波減衰項の係数
!  * スポンジ層の設定(境界付近で波の反射を抑え吸収するための層)
! 
!== Error Handling
!
!== Bugs
!
!== Note
!
!  * この関数は, 基本場が零な変数(速度, エクスナー関数の擾乱)に
!    適用することを想定している. 
!  * 各格子点に対する関数を定義する必要がある
!
!== Future Plans
!
!

module Damping
  !
  !減衰率とその計算を行うためのパッケージ型モジュール
  !  * 音波減衰項の係数
  !  * スポンジ層の設定(境界付近で波の反射を抑え吸収するための層)
  ! 

  !モジュール読み込み
  use gridset, only: DimXMin,       &! x 方向の配列の下限
    &                DimXMax,       &! x 方向の配列の上限
    &                DimZMin,       &! z 方向の配列の下限
    &                DimZMax,       &! z 方向の配列の上限
    &                DelX,          &! x 方向の格子点間隔
    &                DelZ,          &! z 方向の格子点間隔
    &                s_X,           &!X 座標軸(スカラー格子点)
    &                s_Z,           &!Z 座標軸(スカラー格子点)
    &                f_X,           &!X 座標軸(フラックス格子点)
    &                f_Z,           &!Z 座標軸(フラックス格子点)
    &                XMax,          &!X 座標の最大値
    &                ZMax            !Z 座標の最大値
  use timeset,   only: DelTimeShort  !短い時間ステップ
  
  !暗黙の型宣言禁止
  implicit none

  !private 属性を指定
  private 
  
  !関数には public 属性を指定
  public Damping_Init
  public DampSound_Init
  public DampSponge_Init
  public xz_DampSponge_xz
  public xr_DampSponge_xr
  public pz_DampSponge_pz
  public DampSound

  !変数定義
  real(8)  :: DampSound  = 0.0d0   !音波減衰項の減衰係数

  real(8)  :: EFTime     = 100.0d0 !スポンジ層の e-folding time
  real(8)  :: DampDepthH = 0.0d0   !スポンジ層の厚さ(水平方向)
  real(8)  :: DampDepthV = 0.0d0   !スポンジ層の厚さ(鉛直方向)
  real(8), allocatable    :: xz_DampRateH(:,:) !ss 格子点減衰係数(水平方向)
  real(8), allocatable    :: xz_DampRateV(:,:) !ss 格子点減衰係数(鉛直方向)
  real(8), allocatable    :: pz_DampRateH(:,:) !fs 格子点減衰係数(水平方向)
  real(8), allocatable    :: pz_DampRateV(:,:) !fs 格子点減衰係数(鉛直方向)
  real(8), allocatable    :: xr_DampRateH(:,:) !sf 格子点減衰係数(水平方向)
  real(8), allocatable    :: xr_DampRateV(:,:) !sf 格子点減衰係数(鉛直方向)
  
  !値を save する
  save DampSound, EFTime, DampDepthH, DampDepthV
  save xz_DampRateH, xz_DampRateV
  save xr_DampRateH, xr_DampRateV
  save pz_DampRateH, pz_DampRateV
  
contains 
  
!!!------------------------------------------------------------------------!!!
  subroutine Damping_Init( cfgfile ) 
    !
    ! 音波減衰項とスポンジ層の減衰係数の初期化
    ! 

    !暗黙の型宣言禁止
    implicit none

    !変数定義
    character(*), intent(in) :: cfgfile
    real(8)                  :: Alpha
    real(8)                  :: Time
    real(8)                  :: DepthH
    real(8)                  :: DepthV

    !NAMELIST から取得
    NAMELIST /damping/ Alpha, Time, DepthH, DepthV
    open (10, FILE=cfgfile)
    read(10, NML=damping)
    close(10)

    !初期化
    call DampSound_Init( Alpha ) 
    call DampSponge_Init( Time, DepthH, DepthV )


  end subroutine Damping_Init


!!!------------------------------------------------------------------------!!!
  subroutine DampSound_Init( damp ) 
    !
    ! 音波減衰項の減衰係数の初期化
    ! 

    !暗黙の型宣言禁止
    implicit none

    !変数定義
    real(8), intent(in)  :: damp

    !-------------------------------------------------------------------
    ! 音波減衰項の減衰率   DelX ** 2.0 に依存?
    !-------------------------------------------------------------------
    DampSound = Damp * ( DelX ** 2.0d0 ) / DelTimeShort
    
    !-----------------------------------------------------------------    
    ! 値の確認
    !-----------------------------------------------------------------
    write(*,*) "DampSound_init, DampSound : ", DampSound

  end subroutine DampSound_Init


!!!------------------------------------------------------------------------!!!
  subroutine DampSponge_Init( Time, DepthH, DepthV )
    !
    ! スポンジ層の減衰係数を決める
    !

    !暗黙の型宣言禁止
    implicit none

    !入力変数
    real(8), intent(in)   :: Time     !スポンジ層の e-folding time
    real(8), intent(in)   :: DepthH   !スポンジ層の厚さ(水平方向)
    real(8), intent(in)   :: DepthV   !スポンジ層の厚さ(鉛直方向)
    real(8), parameter    :: Pi =3.1415926535897932385d0   !円周率
    integer               :: i, k

    !初期化
    allocate( &
      & xz_DampRateH(DimXMin:DimXMax, DimZMin:DimZMax), &
      & xz_DampRateV(DimXMin:DimXMax, DimZMin:DimZMax), &
      & pz_DampRateH(DimXMin:DimXMax, DimZMin:DimZMax), &
      & pz_DampRateV(DimXMin:DimXMax, DimZMin:DimZMax), &
      & xr_DampRateH(DimXMin:DimXMax, DimZMin:DimZMax), &
      & xr_DampRateV(DimXMin:DimXMax, DimZMin:DimZMax)    )
    xz_DampRateH = 0.0d0
    xz_DampRateV = 0.0d0
    pz_DampRateH = 0.0d0
    pz_DampRateV = 0.0d0
    xr_DampRateH = 0.0d0
    xr_DampRateV = 0.0d0

    !値の入力
    EFTime     = Time
    DampDepthH = DepthH
    DampDepthV = DepthV
    
    !-----------------------------------------------------------------    
    ! スポンジ層の減衰率
    !-----------------------------------------------------------------
    !水平方向の東側・西側境界
    if ( DampDepthH < DelX ) then 
!      write(*,*) "DampDepthH is too thin. DelX is ", DelX   !警告のみ出力

    else
      do i = DimXMin, DimXMax
        !スカラー格子点の西側境界
        if ( s_X(i) < DampDepthH) then 
          xz_DampRateH(i,:) = ((1.0d0 - s_X(i) / DampDepthH) ** 3.0d0) / EFTime
        end if
        
        !フラックス格子点の西側境界
        if ( f_X(i) < DampDepthH) then 
          pz_DampRateH(i,:) = ((1.0d0 - f_X(i) / DampDepthH) ** 3.0d0) / EFTime
         end if
        
        !スカラー格子点の東側境界    
        if ( s_X(i) > ( XMax - DampDepthH ) ) then 
          xz_DampRateH(i,:) = &
            & ((1.0d0 - (XMax - s_X(i)) / DampDepthH) ** 3.0d0) / EFTime 
        end if
        
        !フラックス格子点の東側境界    
        if ( f_X(i) > ( XMax - DampDepthH ) ) then 
          pz_DampRateH(i,:) = &
            & ((1.0d0 - (XMax - f_X(i)) / DampDepthH) ** 3.0d0) / EFTime 
        end if
      end do
    end if
    !sf と ss は X 方向に関しては同じ
    xr_DampRateH  = xz_DampRateH
    
    !鉛直方向の上部境界    
    if ( DampDepthV < DelZ ) then 
!      write(*,*) "DampDepthV is too thin. DelZ is ", DelZ   !警告のみ表示
      
    else
      do k = DimZMin, DimZMax
        !スカラー格子点
        if ( s_Z(k) >= ( ZMax - DampDepthV ) ) then 
          xz_DampRateV(:,k) =  &
            & (1.0d0 - dcos(Pi * (s_Z(k) - ZMax + DampDepthV) / DampDepthV)) &
            &  / EFTime 
        end if
        
        !フラックス格子点
        if ( f_Z(k) >= ( ZMax - DampDepthV ) ) then 
          xr_DampRateV(:,k) =  &
            & (1.0d0 - dcos(Pi * (f_Z(k) - ZMax + DampDepthV)/ DampDepthV)) &
            &  / EFTime 
        end if
      end do
    end if
    !fs と ss は Z 方向に関しては同じ
    pz_DampRateV  = xz_DampRateV
    
    !-----------------------------------------------------------------    
    ! 値の確認
    !-----------------------------------------------------------------
    write(*,*) "DampSponge_Init, EFTime :    ", EFTime 
    write(*,*) "DampSponge_Init, DampDepthH: ", DampDepthH 
    write(*,*) "DampSponge_Init, DampDepthV: ", DampDepthV    

  end subroutine DampSponge_Init
  

!!!------------------------------------------------------------------------!!!
  function xz_DampSponge_xz(xz_Var)
    !
    ! ss 格子点に対するスポンジ層
    !

    !暗黙の型宣言禁止
    implicit none

    !変数定義
    real(8), intent(in)   :: xz_Var(DimXMin:DimXMax, DimZMin:DimZMax)
    real(8)               :: xz_DampSponge_xz(DimXMin:DimXMax, DimZMin:DimZMax)
    
    !スポンジ層によるダンピングを計算
    xz_DampSponge_xz =  ( xz_DampRateH + xz_DampRateV ) * xz_Var
    
  end function xz_DampSponge_xz


!!!------------------------------------------------------------------------!!!
  function xr_DampSponge_xr(xr_Var)
    !
    ! sf 格子点に対するスポンジ層
    !
    
    !暗黙の型宣言禁止
    implicit none

    !変数定義
    real(8), intent(in)   :: xr_Var(DimXMin:DimXMax, DimZMin:DimZMax)
    real(8)               :: xr_DampSponge_xr(DimXMin:DimXMax, DimZMin:DimZMax)

    !スポンジ層によるダンピングを計算  
    xr_DampSponge_xr =  ( xr_DampRateH + xr_DampRateV )* xr_Var
    
  end function xr_DampSponge_xr


!!!------------------------------------------------------------------------!!!
  function pz_DampSponge_pz(pz_Var)
    !
    ! fs 格子点に対するスポンジ層
    !

    !暗黙の型宣言禁止
    implicit none

    !変数定義
    real(8), intent(in)   :: pz_Var(DimXMin:DimXMax, DimZMin:DimZMax)
    real(8)               :: pz_DampSponge_pz(DimXMin:DimXMax, DimZMin:DimZMax)
    
    !スポンジ層によるダンピングを計算  
    pz_DampSponge_pz = ( pz_DampRateH + pz_DampRateV ) * pz_Var
    
  end function pz_DampSponge_pz
  
end module Damping
