Class rad_15m_NLTE
In: radiation/rad_15m_NLTE.f90

(火星大気向け) Non-LTE 放射モデル

Non-NLTE radiation model (for the Mars’ atmosphere)

Note that Japanese and English are described in parallel.

(火星大気向け) Non-LTE 放射モデル

Non-NLTE radiation model (for the Mars’ atmosphere)

References

See Takahashi et al. (2003) for reference

Procedures List

!$ ! RadiationFluxDennouAGCM :放射フラックスの計算
!$ ! RadiationDTempDt :放射フラックスによる温度変化の計算
!$ ! RadiationFluxOutput :放射フラックスの出力
!$ ! RadiationFinalize :終了処理 (モジュール内部の変数の割り付け解除)
!$ ! ———— :————
!$ ! RadiationFluxDennouAGCM :Calculate radiation flux
!$ ! RadiationDTempDt :Calculate temperature tendency with radiation flux
!$ ! RadiationFluxOutput :Output radiation fluxes
!$ ! RadiationFinalize :Termination (deallocate variables in this module)

NAMELIST

!$ ! NAMELIST#rad_15m_NLTE_nml

Methods

Included Modules

dc_types dc_message gridset constants gtool_historyauto timeset constants0 axesset

Public Instance methods

Subroutine :

[Source]

  subroutine Rad15mNLTEInit

    ! ファイル入出力補助
    ! File I/O support
    !
!!$    use dc_iounit, only: FileOpen

    ! NAMELIST ファイル入力に関するユーティリティ
    ! Utilities for NAMELIST file input
    !
!!$    use namelist_util, only: namelist_filename, NmlutilMsg, NmlutilAryValid

    ! ヒストリデータ出力
    ! History data output
    !
    use gtool_historyauto, only: HistoryAutoAddVariable

    ! 座標データ設定
    ! Axes data settings
    !
    use axesset, only: AxnameX, AxnameY, AxnameZ, AxnameR, AxnameT

    ! 宣言文 ; Declaration statements
    !

!!$    integer:: unit_nml        ! NAMELIST ファイルオープン用装置番号.
!!$                              ! Unit number for NAMELIST file open
!!$    integer:: iostat_nml      ! NAMELIST 読み込み時の IOSTAT.
!!$                              ! IOSTAT of NAMELIST read

    ! NAMELIST 変数群
    ! NAMELIST group name
    !
!!$    namelist /rad_15m_NLTE_nml/ &
!!$      & SolarConst
          !
          ! デフォルト値については初期化手続 "rad_15m_NLTE#Rad15mNLTEInit"
          ! のソースコードを参照のこと.
          !
          ! Refer to source codes in the initialization procedure
          ! "rad_15m_NLTE#Rad15mNLTEInit" for the default values.
          !

    if ( rad_15m_NLTE_inited ) return

    ! デフォルト値の設定
    ! Default values settings
    !
!!$    SolarConst      = 1380.0_DP / 1.52_DP**2

    ! NAMELIST の読み込み
    ! NAMELIST is input
    !
!!$    if ( trim(namelist_filename) /= '' ) then
!!$      call FileOpen( unit_nml, &          ! (out)
!!$        & namelist_filename, mode = 'r' ) ! (in)
!!$
!!$      rewind( unit_nml )
!!$      read( unit_nml,                     & ! (in)
!!$        & nml = rad_15m_NLTE_nml,         & ! (out)
!!$        & iostat = iostat_nml )             ! (out)
!!$      close( unit_nml )
!!$
!!$      call NmlutilMsg( iostat_nml, module_name ) ! (in)
!!$    end if


    ! Initialization of modules used in this module
    !


    ! ヒストリデータ出力のためのへの変数登録
    ! Register of variables for history data output
    !
    call HistoryAutoAddVariable( 'DTempDt15mNLTE', (/ AxnameX, AxnameY, AxnameZ, AxnameT /), 'radiative calculation by NLTE process at 15 micron meter', 'K s-1' )
    call HistoryAutoAddVariable( 'DTempDtRadLMerged', (/ AxnameX, AxnameY, AxnameZ, AxnameT /), 'radiative calculation in long wave merged with NLTE heating rate', 'K s-1' )

    ! 印字 ; Print
    !
    call MessageNotify( 'M', module_name, '----- Initialization Messages -----' )
!!$    call MessageNotify( 'M', module_name, 'SolarConst = %f', d = (/ SolarConst /) )
    call MessageNotify( 'M', module_name, '-- version = %c', c1 = trim(version) )

    rad_15m_NLTE_inited = .true.

  end subroutine Rad15mNLTEInit
Subroutine :
xyz_Press(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(in )
xyz_Temp(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(in )
xyz_VirTemp(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(in )
xyz_DTempDt(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(inout)

[Source]

  subroutine rad15mNLTEMergeHR( xyz_Press, xyz_Temp, xyz_VirTemp, xyz_DTempDt )

    ! ヒストリデータ出力
    ! History data output
    !
    use gtool_historyauto, only: HistoryAutoPut

    ! 時刻管理
    ! Time control
    !
    use timeset, only: TimeN                 ! ステップ $ t $ の時刻.
                              ! Time of step $ t $.

    real(DP), intent(in   ) :: xyz_Press  (0:imax-1, 1:jmax, 1:kmax)
    real(DP), intent(in   ) :: xyz_Temp   (0:imax-1, 1:jmax, 1:kmax)
    real(DP), intent(in   ) :: xyz_VirTemp(0:imax-1, 1:jmax, 1:kmax)
    real(DP), intent(inout) :: xyz_DTempDt(0:imax-1, 1:jmax, 1:kmax)


    ! Local variables
    !
    real(DP) :: xyz_Weight        (0:imax-1, 1:jmax, 1:kmax)
    real(DP) :: xyz_DTempDt15mNLTE(0:imax-1, 1:jmax, 1:kmax)


    ! 初期化
    ! Initialization
    !
    if ( .not. rad_15m_NLTE_inited ) then
      call MessageNotify( 'E', module_name, 'This module has not been initialized.' )
    end if


    call rad15mNLTECalcWeight( xyz_Press, xyz_Weight )
!!$    xyz_Weight = 1.0_DP

    call rad15mNLTE( xyz_Press, xyz_Temp, xyz_VirTemp, xyz_DTempDt15mNLTE )

    xyz_DTempDt = xyz_Weight             * xyz_DTempDt + ( 1.0d0 - xyz_Weight ) * xyz_DTempDt15mNLTE


    call HistoryAutoPut( TimeN, 'DTempDt15mNLTE'   , xyz_DTempDt15mNLTE )
    call HistoryAutoPut( TimeN, 'DTempDtRadLMerged', xyz_DTempDt )


  end subroutine rad15mNLTEMergeHR
rad_15m_NLTE_inited
Variable :
rad_15m_NLTE_inited = .false. :logical, save, public
: 初期設定フラグ. Initialization flag

Private Instance methods

a_Table15mNLTEfa
Variable :
a_Table15mNLTEfa( nTable15mNLTE ) :real(DP)
: E(tau)
a_Table15mNLTEsn
Variable :
a_Table15mNLTEsn( nTable15mNLTE ) :real(DP)
: sigma N (non dimension)
module_name
Constant :
module_name = ‘rad_15m_NLTE :character(*), parameter
: モジュールの名称. Module name
nTable15mNLTE
Constant :
nTable15mNLTE = 70 :integer , parameter
Subroutine :
xyz_Press(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(in )
xyz_Temp(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(in )
xyz_VirTemp(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(in )
xyz_DTempDt15mNLTE(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(out)

[Source]

  subroutine rad15mNLTE( xyz_Press, xyz_Temp, xyz_VirTemp, xyz_DTempDt15mNLTE )

    ! 物理定数設定
    ! Physical constants settings
    !
    use constants, only: Grav, CpDry, GasRDry
                              ! $ R $ [J kg-1 K-1].
                              ! 乾燥大気の気体定数.
                              ! Gas constant of air

    real(DP), intent(in ) :: xyz_Press   (0:imax-1, 1:jmax, 1:kmax)
    real(DP), intent(in ) :: xyz_Temp    (0:imax-1, 1:jmax, 1:kmax)
    real(DP), intent(in ) :: xyz_VirTemp (0:imax-1, 1:jmax, 1:kmax)
    real(DP), intent(out) :: xyz_DTempDt15mNLTE(0:imax-1, 1:jmax, 1:kmax)


    ! Local variables
    !
    real(DP) :: xyz_Rho(0:imax-1, 1:jmax, 1:kmax)
    integer  :: i
    integer  :: j
    integer  :: k
    integer  :: l

    real(DP) :: tau
    real(DP) :: e1, e2
    real(DP) :: ltau
    real(DP) :: ramda
    real(DP), parameter :: a10 = 2.44d0
    real(DP), parameter :: g10 = 2.0d0
    real(DP) :: kco2, ko
    real(DP) :: co2nd, ond
    real(DP) :: tmpfac
    real(DP) :: NLTECR
!!$    real(DP) :: xyz_Weight(0:imax-1, 1:jmax, 1:kmax)


    ! Variables for Reference Pressure
    !
    real(DP) :: p0, t0, ond0, f0

    ! kmin     : maximum index used for calculation of
    !          : atmospheric radiative cooling rate
!!$    integer  :: xy_kMin(0:imax-1, 1:jmax)

    real(DP), parameter :: AMU = 1.6605655d-27


    ! 初期化
    ! Initialization
    !
    if ( .not. rad_15m_NLTE_inited ) then
      call MessageNotify( 'E', module_name, 'This module has not been initialized.' )
    end if

!!$    call calcimin(im,iup,jm,km,press,imin)
!!$    xy_kMin = kmax

!!$    call rad15mNLTECalcWeight( &
!!$      & xyz_Press,                   &
!!$!      &  xy_kMin,                    &
!!$      & xyz_Weight                   &
!!$      & )
!!$    xyz_Weight = 1.0_DP

    ! Set Variables for Reference Pressure
    !
    p0 = 1.2d-3 * 1.0d-6 * 1.0d5
    f0 = 0.5d0 * 1.0d-2
!            f0=1.0d0*1.0d-2


    xyz_Rho = xyz_Press / ( GasRDry * xyz_VirTemp )

    do j = 1, jmax
      do i = 0, imax-1

        ! Number Density is in CGS Unit
        !
        if( p0 <= xyz_Press(i,j,kmax) ) then
          t0 = xyz_Temp(i,j,kmax)
          ond0 = xyz_Rho(i,j,kmax) / ( 44.0d0 * AMU ) * ( p0 / xyz_Press(i,j,kmax) ) * f0 * 1.0d-6
        else if( p0 <= xyz_Press(i,j,1) ) then
          search_p : do l = kmax-1, 2, -1
            if( p0 <= xyz_Press(i,j,l) ) exit search_p
          end do search_p
          t0 = ( xyz_Temp(i,j,l+1) - xyz_Temp(i,j,l) ) / log( xyz_Press(i,j,l+1) / xyz_Press(i,j,l) ) * log( p0 / xyz_Press(i,j,l) ) + xyz_Temp(i,j,l)
          ond0 = xyz_Temp(i,j,l) / t0 * xyz_Rho(i,j,l) / ( 44.0d0 * AMU ) * ( p0 / xyz_Press(i,j,l) ) * f0 * 1.0d-6
        else
          write( 6, * ) 'Reference pressure or pressure is inappropriate.'
          write( 6, * ) 'Unexpected Error'
          write( 6, * ) 'Stop'
          stop
        endif

        do k = 1, kmax
          ! cgs unit
          co2nd = xyz_Rho(i,j,k) / ( 44.0d0 * AMU ) * 1.0d-6

!               kco2=1.0d-15
          ! from Lunt et al., 1985
          if( xyz_Temp(i,j,k) < 175.0_DP ) then
            kco2 = 4.2d-12 * exp( -2988.0d0 /175.0d0 + 303930.0d0 / ( 175.0d0 * 175.0d0 ) )
          else
            kco2 = 4.2d-12 * exp( -2988.0d0 / xyz_Temp(i,j,k) + 303930.0d0 / ( xyz_Temp(i,j,k) * xyz_Temp(i,j,k) ) )
          endif

          ! ond0 has already been set in cgs unit
          ond = t0 / xyz_Temp(i,j,k) * ond0 * ( ( xyz_Press(i,j,k) /p0 )**(16.0d0/44.0d0) )

!               ko=1.5d-11*dexp(-800.0d0/ntemp(i,j,k))
          ! from Lopez-Valvelde and Lopez-Puertas, 1994, 
          ! Bougher et al., 1994
          ko = 3.0d-12
          tau = 6.4d-15 * xyz_Press(i,j,k) / ( 44.0d0 * AMU * Grav ) * 1.0d-4


!!$  integer , parameter :: nTable15mNLTE = 70
!!$  real(DP),           :: a_Table15mNLTEsn( nTable15mNLTE )
!!$                                        ! sigma N (non dimension)
!!$  real(DP),           :: a_Table15mNLTEfa( nTable15mNLTE )
!!$                                        ! E(tau)

!----------------------
          search_sn_1 : do l = 2, nTable15mNLTE-1
            if ( tau < a_Table15mNLTEsn(l) ) exit search_sn_1
          end do search_sn_1
          e1 =   ( a_Table15mNLTEfa(l) - a_Table15mNLTEfa(l-1) ) / ( a_Table15mNLTEsn(l) - a_Table15mNLTEsn(l-1) ) * ( tau - a_Table15mNLTEsn(l-1) ) + a_Table15mNLTEfa(l-1)
          if ( e1 > 0.5d0 ) e1 = 0.5d0
          if ( e1 < 0.0d0 ) e1 = 0.0d0
!----------------------
          search_sn_2 : do l = 2, nTable15mNLTE-1
            if ( ( tau / 2.0d0 ) <  a_Table15mNLTEsn(l) ) exit search_sn_2
          end do search_sn_2
          e2 =   ( a_Table15mNLTEfa(l) - a_Table15mNLTEfa(l-1) ) / ( a_Table15mNLTEsn(l) - a_Table15mNLTEsn(l-1) ) * ( ( tau / 2.0d0 ) - a_Table15mNLTEsn(l-1) ) + a_Table15mNLTEfa(l-1)
          if ( e2 > 0.5d0 ) e2 = 0.5d0
          if ( e2 < 0.0d0 ) e2 = 0.0d0
!----------------------
          ltau = e1 + e2
          ramda = a10 / ( a10 + kco2 * co2nd + ko * ond )
          tmpfac = 0.5d0 * ramda * ltau / ( 1.0d0 - ramda + 0.5d0 * ramda * ltau )

          NLTECR = 1.33d-13 * g10 * exp( -960.0d0 / xyz_Temp(i,j,k) ) * co2nd * ( kco2 * co2nd + ko * ond ) * tmpfac * 1.0d-1 / ( xyz_Rho(i,j,k) * CpDry )

!!$          xyz_DTempDt15mNLTE(i,j,k) = &
!!$            & ( 1.0d0 - xyz_Weight(i,j,k) ) * ( -NLTECR )
          xyz_DTempDt15mNLTE(i,j,k) = - NLTECR
        end do

      end do
    end do

  end subroutine rad15mNLTE
Subroutine :
xyz_Press(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(in )
:
!$ integer , intent(in ) :xy_kMin (0:imax-1, 1:jmax)
xyz_Weight(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(out)

[Source]

  subroutine rad15mNLTECalcWeight( xyz_Press, xyz_Weight )

    ! 物理・数学定数設定
    ! Physical and mathematical constants settings
    !
    use constants0, only: PI                    ! $ \pi $ .
                              ! 円周率.  Circular constant

    real(DP), intent(in ) :: xyz_Press (0:imax-1, 1:jmax, 1:kmax)
!!$    integer , intent(in ) :: xy_kMin   (0:imax-1, 1:jmax)
    real(DP), intent(out) :: xyz_Weight(0:imax-1, 1:jmax, 1:kmax)


    ! Local variables
    !
    integer :: i
    integer :: j
    integer :: k

    ! 初期化
    ! Initialization
    !
    if ( .not. rad_15m_NLTE_inited ) then
      call MessageNotify( 'E', module_name, 'This module has not been initialized.' )
    end if

    do k = 1, kmax
      do j = 1, jmax
        do i = 0, imax-1

!               weight(i,j,k)=(atan(2.0d0 &
!                    *dlog(dsqrt(press(i,j,k)*press(i+1,j,k)) &
!                    /dsqrt(press(imin+4,j,k)*press(imin+4+1,j,k)))) &
!                    *1.2d0 &
!                    +pi/2.0d0)/pi

          xyz_Weight(i,j,k) = ( atan( 2.0d0 * log( xyz_Press(i,j,k) / ( 1.0d-2 * exp( 2.0d0 ) ) ) ) * 1.2d0 + Pi / 2.0d0 ) / Pi
          xyz_Weight(i,j,k) = max( xyz_Weight(i,j,k), 0.0d0 )
          xyz_Weight(i,j,k) = min( xyz_Weight(i,j,k), 1.0d0 )
        end do
      end do
    end do


  end subroutine rad15mNLTECalcWeight
Subroutine :
xyz_Press(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(in )
xy_kMin(0:imax-1, 1:jmax) :integer , intent(out)

[Source]

  subroutine rad15mNLTECalckMin( xyz_Press, xy_kMin )

    real(DP), intent(in ) :: xyz_Press(0:imax-1, 1:jmax, 1:kmax)
    integer , intent(out) :: xy_kMin  (0:imax-1, 1:jmax)

    integer :: i
    integer :: j
    integer :: k

    ! 初期化
    ! Initialization
    !
    if ( .not. rad_15m_NLTE_inited ) then
      call MessageNotify( 'E', module_name, 'This module has not been initialized.' )
    end if

    ! set minimum index used for calculation of
    ! atmospheric radiative cooling rate
!    imin=im-25

    do j = 1, jmax
      do i = 0, imax-1
        find_kmin : do k = kmax, 1, -1
          if ( xyz_Press(i,j,k) > 1.0d-2 ) exit find_kmin
        end do find_kmin
        xy_kMin(i,j) = k
      end do
    end do

  end subroutine rad15mNLTECalckMin
version
Constant :
version = ’$Name: dcpam5-20130921 $’ // ’$Id: rad_15m_NLTE.f90,v 1.1 2012-11-15 03:30:10 yot Exp $’ :character(*), parameter
: モジュールのバージョン Module version