!---------------------------------------------------------------------
!     Copyright (C) GFD Dennou Club, 2006. All rights reserved.
!---------------------------------------------------------------------
!= Module Radiation
!
!   * Developer: SUGIYAMA Ko-ichiro
!   * Version: $Id: $
!   * Tag Name: $Name:  $
!   * Change History: 
!
!== Overview
!
!モデルの放射過程を計算するためのパッケージ型モジュール
!具体的には以下の項を計算するための関数を格納する.  
!  * 一様冷却
!
!== Error Handling
!
!== Bugs
!
!== Note
!
!
!== Future Plans
!
!

module Radiation
  !
  !モデルの放射過程を計算するためのパッケージ型モジュール
  !具体的には以下の項を計算するための関数を格納する.  
  !  * 一様冷却
  !

  !モジュール読み込み
  use GridSet,  only: DimXMin,     & !x 方向の配列の下限
    &                 DimXMax,     & !x 方向の配列の上限
    &                 DimZMin,     & !z 方向の配列の下限
    &                 DimZMax,     & !z 方向の配列の上限
    &                 s_Z            !Z 座標軸(スカラー格子点)
  use TimeSet,  only: DayTime        ! 1 日の長さ [s]
  use basicset, only: xz_ExnerBasicZ !エクスナー関数の基本場
  
  !暗黙の型宣言禁止
  implicit none

  !private 属性
  private

  !関数を public にする. 
  public xz_RadCoolConst, Radiation_init

  !変数定義
  real(8) :: RadHeatRate = 0.0d0   !一様放射加熱率 [K/day]
  real(8) :: RadHeight   = 0.0d0   !放射加熱が存在する領域の最大高度 [m]

  !save 属性
  save RadHeatRate, RadHeight

contains

!!!------------------------------------------------------------------------!!!
  subroutine Radiation_init(cfgfile)
    !
    !NAMELIST から放射強制の設定を取得
    !

    !暗黙の型宣言禁止
    implicit none

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

    ! NAMELIST から情報を取得
    NAMELIST /radiation/ RadHeatRate, RadHeight

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

    write(*,*) "Radiation, RadHeatRate: ", RadHeatRate
    write(*,*) "Radiation, RadHeight:   ", RadHeight
    
  end subroutine Radiation_init


!!!------------------------------------------------------------------------!!!
  function xz_RadCoolConst(xz_Exner)
    !
    ! 温位の放射強制項. 
    ! 地表面から RadHeight で指定された高度までの間で一様放射冷却を与える. 
    !
    
    !暗黙の型宣言を禁止
    implicit none

    !変数定義
    real(8), intent(in)   :: xz_Exner(DimXMin:DimXMax, DimZMin:DimZMax)
                                       !エクスナー関数の擾乱場
    real(8)               :: xz_RadCoolConst(DimXMin:DimXMax, DimZMin:DimZMax)
                                       !放射強制
    integer               :: k         !ループ変数

    xz_RadCoolConst = 0.0d0
    do k = DimZMin, DimZMax
      if ( s_Z(k) <= RadHeight ) then
        xz_RadCoolConst(DimXMin:DimXMax, k) = RadHeatRate   &
          & / (     &
          &     (     &
          &         xz_ExnerBasicZ(DimXMin:DimXMax, k)   &
          &       + xz_Exner(DimXMin:DimXMax, k)      &
          &     ) * DayTime    &
          &   )
      else
        xz_RadCoolConst(DimXMin:DimXMax, k) = - 0.0d0
      end if
    end do
    
  end function xz_RadCoolConst
  
end module Radiation
