!---------------------------------------------------------------------
!     Copyright (C) GFD Dennou Club, 2004, 2006. All rights reserved.
!---------------------------------------------------------------------
!= Module GridSet
!
!   * Developer: SUGIYAMA Ko-ichiro 
!   * Version: $Id: $ 
!   * Tag Name: $Name:  $
!   * Change History: 
!
!== Overview 
!
!引数に与えられた NAMELIST ファイルから, 格子点情報を取得し, 
!保管するための変数型モジュール
!
!== Error Handling
!
!== Known Bugs
!
!== Note
!
!== Future Plans
!

module gridset
  !
  !引数に与えられた NAMELIST ファイルから, 格子点情報を取得し, 
  !保管するための変数型モジュール
  !

  !暗黙の型宣言禁止
  implicit none

  !save 属性
  save

  !公開変数
  real(8)               :: Xmin, Xmax    ! x 座標の始点・終点
  real(8)               :: Zmin, Zmax    ! z 座標の始点・終点
  real(8)               :: DelX, DelZ    !格子間隔
  integer               :: NX, NZ        !格子点数
  integer, parameter    :: Margin = 5    !境界のグリッド数
  integer               :: SpcNum        !化学種の数
  integer               :: DimXMin       ! x 方向の配列の下限
  integer               :: DimXMax       ! x 方向の配列の上限
  integer               :: DimZMin       ! z 方向の配列の下限
  integer               :: DimZMax       ! z 方向の配列の上限
  integer               :: RegXMin       ! x 方向の物理領域の下限
  integer               :: RegXMax       ! x 方向の物理領域の上限
  integer               :: RegZMin       ! z 方向の物理領域の下限
  integer               :: RegZMax       ! z 方向の物理領域の上限
  integer               :: FileNX        !ファイル出力用
  integer               :: FileNZ        !ファイル出力用
  integer               :: FileXMin      !ファイル出力用
  integer               :: FileXMax      !ファイル出力用
  integer               :: FileZMin      !ファイル出力用
  integer               :: FileZMax      !ファイル出力用
  real(8), allocatable  :: s_X(:)        !X 座標軸(スカラー格子点)
  real(8), allocatable  :: f_X(:)        !X 座標軸(ベクトル格子点)
  real(8), allocatable  :: s_Z(:)        !Z 座標軸(スカラー格子点)
  real(8), allocatable  :: f_Z(:)        !Z 座標軸(ベクトル格子点)

contains

  subroutine gridset_init(cfgfile)
    !
    !NAMELIST から情報を得て, 格子点を計算する
    !

    !モジュール読み込み
    use DebugSet, only: DebugOn
    use chemcalc, only: ChemCalcDim_init !配列の大きさの初期化
    
    !暗黙の型宣言禁止
    implicit none

    !入力変数
    character(*), intent(in) :: cfgfile

    !内部変数
    integer            :: i, k
    integer, parameter :: kind = 8      !精度を表す

    !-----------------------------------------------------------------
    ! NAMELIST から情報を取得
    !-----------------------------------------------------------------
    NAMELIST /gridset/ NX, NZ, Xmin, Xmax, Zmin, Zmax, SpcNum

    open (10, FILE=cfgfile)
    read(10, NML=gridset)
    close(10)

    !-----------------------------------------------------------------
    ! 格子点間隔計算
    !-----------------------------------------------------------------
    DelX = (Xmax - Xmin) / real(NX, kind)
    DelZ = (Zmax - Zmin) / real(NZ, kind)

    !-----------------------------------------------------------------
    ! 物理的に意味のある領域の上限・下限を決める
    !-----------------------------------------------------------------;
    RegXMin = 0
    RegXMax = NX
    RegZMin = 0
    RegZMax = NZ

    !-----------------------------------------------------------------
    ! 配列の上限・下限を決める
    !-----------------------------------------------------------------
    DimXMin = RegXMin - Margin
    DimXMax = RegXMax + Margin
    DimZMin = RegZMin - Margin
    DimZMax = RegZMax + Margin

    !-----------------------------------------------------------------
    ! グリッドの設定
    !-----------------------------------------------------------------
    allocate( f_X(DimXMin:DimXMax), f_Z(DimZMin:DimZMax) )
    do i = DimXMin, DimXMax
      f_X(i) = Xmin + DelX * ( real(i, kind) - RegXMin )
    end do
    do k = DimZMin, DimZMax
      f_Z(k) = Zmin + DelZ * ( real(k, kind) - RegZMin )
    end do
    
    !-----------------------------------------------------------------
    ! 半格子ずれたグリッドの設定. 等間隔を前提にしている.
    !-----------------------------------------------------------------
    allocate( s_X(DimXMin: DimXMax), s_Z(DimZMin: DimZMax) )
    s_X(DimXMin+1:DimXMax) =           &
      &  (                             &
      &     f_X(DimXMin+1:DimXMax)     &
      &   + f_X(DimXMin:DimXMax-1)     &
      &  ) * 5.0d-1
    s_X(DimXMin) = s_X(DimXMin+1) - DelX
    
    s_Z(DimZMin+1:DimZMax) =           &
      &  (                             &
      &     f_Z(DimZMin+1:DimZMax)     &
      &   + f_Z(DimZMin:DimZMax-1)     &
      &  ) * 5.0d-1
    s_Z(DimZMin) = s_Z(DimZMin+1) - DelZ

    !-----------------------------------------------------------------
    ! デバッグモードか否かで, ファイル出力に用いる変数の大きさを変える
    !-----------------------------------------------------------------
    if (DebugOn) then
      FileNX   = size(f_X, 1)
      FileNZ   = size(f_Z, 1)
      FileXMin = DimXMin
      FileXMax = DimXMax
      FileZMin = DimZMin
      FileZMax = DimZMax
    else
      FileNX   = NX
      FileNZ   = NZ
      FileXMin = RegXMin + 1
      FileXMax = RegXMax
      FileZMin = RegZMin + 1
      FileZMax = RegZMax
    end if

    !-----------------------------------------------------------------    
    ! 確認
    !-----------------------------------------------------------------
    write(*,*) "gridset_init, XMin   ", XMin
    write(*,*) "gridset_init: XMax   ", XMax
    write(*,*) "gridset_init: ZMin   ", ZMin
    write(*,*) "gridset_init: ZMax   ", ZMax
    write(*,*) "gridset_init: DelX   ", DelX
    write(*,*) "gridset_init: DelZ   ", DelZ
    write(*,*) "gridset_init: NX     ", NX
    write(*,*) "gridset_init: NZ     ", NZ
    write(*,*) "gridset_init: SpcNum ", SpcNum
    write(*,*) "gridset_init: Margin ", Margin
    write(*,*) "gridset_init: DimXMin", DimXMin
    write(*,*) "gridset_init: DimXMax", DimXMax
    write(*,*) "gridset_init: DimZMin", DimZMin
    write(*,*) "gridset_init: DimZMax", DimZMax
  end subroutine gridset_init
  
end module gridset
