Class rad_simple_LW
In: radiation/rad_simple_LW.f90

放射フラックス (簡単長波バンドモデル)

Radiation flux (simple longwave band model)

Note that Japanese and English are described in parallel.

温度, 比湿, 気圧から, 放射フラックスを計算する放射モデルです.

This is a radiation model that calculates radiation flux from temperature, specific humidity, and air pressure.

References

Procedures List

RadSimpleLWFlux :放射フラックスの計算
!$ ! RadSimpleLWFinalize :終了処理 (モジュール内部の変数の割り付け解除)
——————- :————
RadSimpleLWFlux :Calculate radiation flux
!$ ! RadSimpleLWFinalize :Termination (deallocate variables in this module)

NAMELIST

NAMELIST#rad_Simple_LW_nml

Methods

Included Modules

gridset dc_types namelist_util dc_message gtool_history constants0 planck_func rad_rte_nonscat fileset axesset timeset dc_calendar dc_iounit dc_string

Public Instance methods

Subroutine :
xy_SurfAlbedo(0:imax-1, 1:jmax) :real(DP), intent(in)
: Surface albedo
xy_SurfEmis(0:imax-1, 1:jmax) :real(DP), intent(in)
: 地表面射出率. Surface emissivity
xyr_Press(0:imax-1, 1:jmax, 0:kmax) :real(DP), intent(in)
: $ P $ . 圧力. Pressure
xyz_Press(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(in)
: $ P $ . 圧力. Pressure
xyz_Temp(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(in)
: $ T $ . 温度. Temperature
xy_SurfTemp(0:imax-1, 1:jmax) :real(DP), intent(in)
: 地表面温度. Surface temperature
xyz_DelAtmMass(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(in)
: Atmospheric mass of layers
xyz_QH2OVap(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(in)
: Specific humidity
xyr_RadLUwFlux(0:imax-1, 1:jmax, 0:kmax) :real(DP), intent(out)
: 長波フラックス. Upward longwave flux
xyr_RadLDwFlux(0:imax-1, 1:jmax, 0:kmax) :real(DP), intent(out)
: 長波フラックス. Downward longwave flux
xyra_DelRadLUwFlux(0:imax-1, 1:jmax, 0:kmax, 0:1) :real(DP), intent(out)
: 長波地表温度変化.
xyra_DelRadLDwFlux(0:imax-1, 1:jmax, 0:kmax, 0:1) :real(DP), intent(out)
: 長波地表温度変化.

長波フラックスの計算

Calculate long wave flux

[Source]

  subroutine RadSimpleLWFlux( xy_SurfAlbedo, xy_SurfEmis, xyr_Press, xyz_Press, xyz_Temp, xy_SurfTemp, xyz_DelAtmMass, xyz_QH2OVap, xyr_RadLUwFlux, xyr_RadLDwFlux, xyra_DelRadLUwFlux, xyra_DelRadLDwFlux )
    !
    ! 長波フラックスの計算
    !
    ! Calculate long wave flux
    !

    ! モジュール引用 ; USE statements
    !

    ! 物理・数学定数設定
    ! Physical and mathematical constants settings
    !
    use constants0, only: PI, StB                   ! $ \sigma_{SB} $ . 
                              ! ステファンボルツマン定数. 
                              ! Stefan-Boltzmann constant

    ! プランク関数の計算
    ! Calculate Planck function
    !
    use planck_func, only : Integ_PF_GQ_Array3D, Integ_PF_GQ_Array2D, Integ_DPFDT_GQ_Array2D

    ! 散乱を無視した放射伝達方程式
    ! Radiative transfer equation without considering scattering
    !
    use rad_rte_nonscat, only : RadRTENonScatMonoSemiAnal


    ! 宣言文 ; Declaration statements
    !
    real(DP), intent(in):: xy_SurfAlbedo     (0:imax-1, 1:jmax)
                              !
                              ! Surface albedo
    real(DP), intent(in):: xy_SurfEmis       (0:imax-1, 1:jmax)
                              ! 地表面射出率.
                              ! Surface emissivity
    real(DP), intent(in):: xyr_Press         (0:imax-1, 1:jmax, 0:kmax)
                              ! $ P $ .     圧力. Pressure
    real(DP), intent(in):: xyz_Press         (0:imax-1, 1:jmax, 1:kmax)
                              ! $ P $ .     圧力. Pressure
    real(DP), intent(in):: xyz_Temp          (0:imax-1, 1:jmax, 1:kmax)
                              ! $ T $ .     温度. Temperature
    real(DP), intent(in):: xy_SurfTemp       (0:imax-1, 1:jmax)
                              ! 地表面温度. 
                              ! Surface temperature
    real(DP), intent(in):: xyz_DelAtmMass    (0:imax-1, 1:jmax, 1:kmax)
                              ! 
                              ! Atmospheric mass of layers
    real(DP), intent(in):: xyz_QH2OVap       (0:imax-1, 1:jmax, 1:kmax)
                              !
                              ! Specific humidity
    real(DP), intent(out):: xyr_RadLUwFlux     (0:imax-1, 1:jmax, 0:kmax)
                              ! 長波フラックス. 
                              ! Upward longwave flux
    real(DP), intent(out):: xyr_RadLDwFlux     (0:imax-1, 1:jmax, 0:kmax)
                              ! 長波フラックス. 
                              ! Downward longwave flux
    real(DP), intent(out):: xyra_DelRadLUwFlux (0:imax-1, 1:jmax, 0:kmax, 0:1)
                              ! 長波地表温度変化. 
                              ! 
    real(DP), intent(out):: xyra_DelRadLDwFlux (0:imax-1, 1:jmax, 0:kmax, 0:1)
                              ! 長波地表温度変化. 
                              ! 

    ! 作業変数
    ! Work variables
    !
    real(DP) :: xyr_Temp    (0:imax-1, 1:jmax, 0:kmax)
    real(DP) :: xyr_IntPF   (0:imax-1, 1:jmax, 0:kmax)
    real(DP) :: xy_SurfIntPF(0:imax-1, 1:jmax)
    real(DP) :: xy_IntDPFDT0(0:imax-1, 1:jmax)
    real(DP) :: xy_IntDPFDT1(0:imax-1, 1:jmax)

    real(DP) :: xyz_DelOptDepDryCom(0:imax-1, 1:jmax, 1:kmax)
    real(DP) :: xyz_DelOptDepH2OVap(0:imax-1, 1:jmax, 1:kmax)

    real(DP) :: xyr_RadUwFlux     (0:imax-1, 1:jmax, 0:kmax)
                              ! 長波フラックス. 
                              ! Upward longwave flux
    real(DP) :: xyr_RadDwFlux     (0:imax-1, 1:jmax, 0:kmax)
                              ! 長波フラックス. 
                              ! Downward longwave flux
    real(DP) :: xyra_DelRadUwFlux (0:imax-1, 1:jmax, 0:kmax, 0:1)
                              ! 長波地表温度変化. 
                              ! 
    real(DP) :: xyra_DelRadDwFlux (0:imax-1, 1:jmax, 0:kmax, 0:1)
                              ! 長波地表温度変化. 
                              ! 

    real(DP) :: WNs
    real(DP) :: WNe

    integer:: k               ! 鉛直方向に回る DO ループ用作業変数
                              ! Work variables for DO loop in vertical direction
    integer:: n               ! 波長について回る DO ループ用作業変数
                              ! Work variables for DO loop in wavenumber bands

    ! 実行文 ; Executable statement
    !


    k = 0
    xyr_Temp(:,:,k) = ( xyz_Temp (:,:,k+2) - xyz_Temp (:,:,k+1) ) / log( xyz_Press(:,:,k+2) / xyz_Press(:,:,k+1) ) * log( xyr_Press(:,:,k  ) / xyz_Press(:,:,k+1) ) + xyz_Temp(:,:,k+1)
    do k = 1, kmax-1
      xyr_Temp(:,:,k) = ( xyz_Temp (:,:,k+1) - xyz_Temp (:,:,k) ) / log( xyz_Press(:,:,k+1) / xyz_Press(:,:,k) ) * log( xyr_Press(:,:,k  ) / xyz_Press(:,:,k) ) + xyz_Temp(:,:,k)
    end do
    k = kmax
    xyr_Temp(:,:,k) = xyz_Temp(:,:,k)

    !   Initialization
    !
    xyr_RadLUwFlux     = 0.0_DP
    xyr_RadLDwFlux     = 0.0_DP
    xyra_DelRadLUwFlux = 0.0_DP
    xyra_DelRadLDwFlux = 0.0_DP
    !
    LOOP_BAND_RTE : do n = 1, nbmax


      ! $ \pi B $, $ \pi DBDT $ の計算
      ! Calculate $ \pi B $ and $ \pi DBDT $
      !
      if ( nbmax == 1 ) then

        xy_SurfIntPF = xy_SurfEmis * StB * ( xy_SurfTemp**4 )
        xyr_IntPF    =               StB * ( xyr_Temp**4 )
        xy_IntDPFDT0 = xy_SurfEmis * 4.0_DP * xy_SurfIntPF / xy_SurfTemp
!!$        xy_IntDPFDT1 =               4.0_DP * xyz_IntPF(:,:,1) / xyz_Temp(:,:,1)
        xy_IntDPFDT1 = 0.0_DP

      else

        WNs = a_WNBnds(n-1)
        WNe = a_WNBnds(n  )
        call Integ_PF_GQ_Array3D( WNs, WNe, NumGaussNode, 0, imax-1, 1, jmax, 0, kmax, xyr_Temp, xyr_IntPF )
        call Integ_PF_GQ_Array2D( WNs, WNe, NumGaussNode, 0, imax-1, 1, jmax, xy_SurfTemp, xy_SurfIntPF )
!!$        call Integ_DPFDT_GQ_Array2D(             &
!!$          & WNs, WNe, NumGaussNode,              & ! (in )
!!$          & 0, imax-1, 1, jmax, xyz_Temp(:,:,1), & ! (in )
!!$          & xy_IntDPFDT1                         & ! (out)
!!$          & )
        xy_IntDPFDT1 = 0.0_DP
        call Integ_DPFDT_GQ_Array2D( WNs, WNe, NumGaussNode, 0, imax-1, 1, jmax, xy_SurfTemp, xy_IntDPFDT0 )

        xy_SurfIntPF = xy_SurfEmis * PI * xy_SurfIntPF
        xyr_IntPF    =               PI * xyr_IntPF
        xy_IntDPFDT0 = xy_SurfEmis * PI * xy_IntDPFDT0
        xy_IntDPFDT1 =               PI * xy_IntDPFDT1

      end if


      ! 光学的厚さの計算
      ! Calculate optical depth
      !
      xyz_DelOptDepDryCom = a_AbsCoefDryCom(n) * ( xyz_Press / a_RefPressDryCom(n) )**a_PressScaleIndDryCom(n) * xyz_DelAtmMass * RadActDryComMMR
      xyz_DelOptDepH2OVap = a_AbsCoefH2OVap(n) * ( xyz_Press / a_RefPressH2OVap(n) )**a_PressScaleIndH2OVap(n) * xyz_DelAtmMass * xyz_QH2OVap


      call RadRTENonScatMonoSemiAnal( xy_SurfAlbedo, xyr_IntPF, xy_SurfIntPF, xy_IntDPFDT1, xy_IntDPFDT0, ( xyz_DelOptDepDryCom + xyz_DelOptDepH2OVap ), xyr_RadUwFlux, xyr_RadDwFlux, xyra_DelRadUwFlux, xyra_DelRadDwFlux )


      xyr_RadLUwFlux     = xyr_RadLUwFlux     + xyr_RadUwFlux
      xyr_RadLDwFlux     = xyr_RadLDwFlux     + xyr_RadDwFlux
      xyra_DelRadLUwFlux = xyra_DelRadLUwFlux + xyra_DelRadUwFlux
      xyra_DelRadLDwFlux = xyra_DelRadLDwFlux + xyra_DelRadDwFlux

    end do LOOP_BAND_RTE


  end subroutine RadSimpleLWFlux
Subroutine :

rad_simple_LW モジュールの初期化を行います. NAMELIST#rad_simple_LW_nml の読み込みはこの手続きで行われます.

"rad_simple_LW" module is initialized. "NAMELIST#rad_simple_LW_nml" is loaded in this procedure.

This procedure input/output NAMELIST#rad_simple_LW_nml .

[Source]

  subroutine RadSimpleLWInit
    !
    ! rad_simple_LW モジュールの初期化を行います. 
    ! NAMELIST#rad_simple_LW_nml の読み込みはこの手続きで行われます. 
    !
    ! "rad_simple_LW" module is initialized. 
    ! "NAMELIST#rad_simple_LW_nml" is loaded in this procedure. 
    !

    ! モジュール引用 ; USE statements
    !

    ! 出力ファイルの基本情報
    ! Basic information for output files
    ! 
    use fileset, only: FileTitle, FileSource, FileInstitution
                              ! データファイルを最終的に変更した組織/個人. 
                              ! Institution or person that changes data files for the last time

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

    ! 座標データ設定
    ! Axes data settings
    !
    use axesset, only: x_Lon, x_Lon_Weight, y_Lat, y_Lat_Weight, z_Sigma, r_Sigma, z_DelSigma
                              ! $ \Delta \sigma $ (整数). 
                              ! $ \Delta \sigma $ (Full)

    ! 時刻管理
    ! Time control
    !
    use timeset, only: RestartTime           ! リスタート開始時刻. 
                              ! Retart time of calculation


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

    ! 暦と日時の取り扱い
    ! Calendar and Date handler
    !
    use dc_calendar, only: DCCalConvertByUnit

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

    ! 種別型パラメタ
    ! Kind type parameter
    !
    use dc_types, only: STDOUT ! 標準出力の装置番号. Unit number of standard output

    ! 文字列操作
    ! Character handling
    !
    use dc_string, only: toChar

    ! 散乱を無視した放射伝達方程式
    ! Radiative transfer equation without considering scattering
    !
    use rad_rte_nonscat, only : RadRTENonScatInit

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

    real(DP) :: WNBnds             (0:MaxNmlArySize)
    real(DP) :: AbsCoefDryCom      (1:MaxNmlArySize)
    real(DP) :: PressScaleIndDryCom(1:MaxNmlArySize)
    real(DP) :: RefPressDryCom     (1:MaxNmlArySize)
    real(DP) :: AbsCoefH2OVap      (1:MaxNmlArySize)
    real(DP) :: PressScaleIndH2OVap(1:MaxNmlArySize)
    real(DP) :: RefPressH2OVap     (1:MaxNmlArySize)

    integer  :: n

    ! NAMELIST 変数群
    ! NAMELIST group name
    !
    namelist /rad_simple_LW_nml/ nbmax, WNBnds, AbsCoefDryCom, PressScaleIndDryCom, RefPressDryCom, AbsCoefH2OVap, PressScaleIndH2OVap, RefPressH2OVap, RadActDryComMMR, DiffFact, NumGaussNode
          !
          ! デフォルト値については初期化手続 "rad_DennouAGCM#RadInit" 
          ! のソースコードを参照のこと. 
          !
          ! Refer to source codes in the initialization procedure
          ! "rad_DennouAGCM#RadInit" for the default values. 
          !

    ! 実行文 ; Executable statement
    !

    if ( rad_simple_LW_inited ) return


    ! デフォルト値の設定
    ! Default values settings
    !

    ! 長波フラックス用情報
    ! Information for long wave flux
    !

    nbmax                  = 1

    WNBnds                 = -999.9_DP
    AbsCoefDryCom          = -999.9_DP
    PressScaleIndDryCom    = -999.9_DP
    RefPressDryCom         = -999.9_DP
    AbsCoefH2OVap          = -999.9_DP
    PressScaleIndH2OVap    = -999.9_DP
    RefPressH2OVap         = -999.9_DP

    AbsCoefDryCom      (1:nbmax) = (/ 5.0d-5 /)
    PressScaleIndDryCom(1:nbmax) = (/ 0.0d0  /)
    RefPressDryCom     (1:nbmax) = (/ 1.0d5  /)

    AbsCoefH2OVap      (1:nbmax) = (/ 1.0d-2 /)
    PressScaleIndH2OVap(1:nbmax) = (/ 0.0d0  /)
    RefPressH2OVap     (1:nbmax) = (/ 1.0d5  /)

    RadActDryComMMR       = 1.0_DP

    DiffFact        = 1.66_DP

    NumGaussNode    = 5


    ! NAMELIST の読み込み
    ! NAMELIST is input
    !
    if ( trim(namelist_filename) /= '' ) then
      call FileOpen( unit_nml, namelist_filename, mode = 'r' ) ! (in)

      rewind( unit_nml )
      read( unit_nml, nml = rad_simple_LW_nml, iostat = iostat_nml )            ! (out)
      close( unit_nml )

      call NmlutilMsg( iostat_nml, module_name ) ! (in)
    end if

    if ( nbmax > MaxNmlArySize ) then
      call MessageNotify( 'E', module_name, 'nbmax = %d > %d', i = (/ nbmax, MaxNmlArySize /) )
    end if

    a_WNBnds              = WNBnds * 100.0_DP    ! Convert from cm-1 to m-1
!!$    a_AbsCoefDryCom = 5.0d-5
    a_AbsCoefDryCom       = AbsCoefDryCom
    a_PressScaleIndDryCom = PressScaleIndDryCom
    a_RefPressDryCom      = RefPressDryCom

!!$    a_AbsCoefH2OVap = 1.0d-2
    a_AbsCoefH2OVap = AbsCoefH2OVap
    a_PressScaleIndH2OVap = PressScaleIndH2OVap
    a_RefPressH2OVap      = RefPressH2OVap


    ! Initialization of modules used in this module
    !

    ! 散乱を無視した放射伝達方程式
    ! Radiative transfer equation without considering scattering
    !
    call RadRTENonScatInit


    ! 印字 ; Print
    !
    call MessageNotify( 'M', module_name, '----- Initialization Messages -----' )

!!$    call MessageNotify( 'M', module_name, 'DelTime:' )
!!$    call MessageNotify( 'M', module_name, '  DelTime  = %f [%c]', &
!!$      & d = (/ DelTimeValue /), c1 = trim( DelTimeUnit ) )

    call MessageNotify( 'M', module_name, 'nbmax         = %d', i = (/ nbmax /) )

!!$    call MessageNotify( 'M', module_name, 'WNBnds              = (/ %*r /)', &
!!$      & r = real( a_WNBnds(0:nbmax) ), n = (/ nbmax+1 /) )
!!$    call MessageNotify( 'M', module_name, 'AbsCoefDryCom       = (/ %*r /)', &
!!$      & r = real( a_AbsCoefDryCom(1:nbmax) ), n = (/ nbmax /) )
!!$    call MessageNotify( 'M', module_name, 'PressScaleIndDryCom = (/ %*r /)', &
!!$      & r = real( a_PressScaleIndDryCom(1:nbmax) ), n = (/ nbmax /) )
!!$    call MessageNotify( 'M', module_name, 'RefPressDryCom      = (/ %*r /)', &
!!$      & r = real( a_RefPressDryCom(1:nbmax) ), n = (/ nbmax /) )
!!$    call MessageNotify( 'M', module_name, 'AbsCoefH2OVap       = (/ %*r /)', &
!!$      & r = real( a_AbsCoefH2OVap(1:nbmax) ), n = (/ nbmax /) )
!!$    call MessageNotify( 'M', module_name, 'PressScaleIndH2OVap = (/ %*r /)', &
!!$      & r = real( a_PressScaleIndH2OVap(1:nbmax) ), n = (/ nbmax /) )
!!$    call MessageNotify( 'M', module_name, 'RefPressH2OVap      = (/ %*r /)', &
!!$      & r = real( a_RefPressH2OVap(1:nbmax) ), n = (/ nbmax /) )

    do n = 1, nbmax
      call MessageNotify( 'M', module_name, '  %d : %f %f %f %f %f %f %f %f', i = (/ n /), d = (/ a_WNBnds(n-1)*1.0e-2, a_WNBnds(n)*1.0e-2, a_AbsCoefDryCom(n), a_PressScaleIndDryCom(n), a_RefPressDryCom(n), a_AbsCoefH2OVap(n), a_PressScaleIndH2OVap(n), a_RefPressH2OVap(n) /) )
    end do

    call MessageNotify( 'M', module_name, 'RadActDryComMMR = %f', d = (/ RadActDryComMMR /) )
    call MessageNotify( 'M', module_name, 'DiffFact        = %f', d = (/ DiffFact /) )
!
    call MessageNotify( 'M', module_name, '-- version = %c', c1 = trim(version) )

    rad_simple_LW_inited = .true.

  end subroutine RadSimpleLWInit

Private Instance methods

DiffFact
Variable :
DiffFact :real(DP), save
: Diffusivity factor
NumGaussNode
Variable :
NumGaussNode :integer , save
Subroutine :
xy_SurfEmis(0:imax-1, 1:jmax) :real(DP), intent(in)
: 地表面射出率. Surface emissivity
xyz_Temp(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(in)
: $ T $ . 温度. Temperature
xyz_Press(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(in)
: $ T $ . 温度. Temperature
xy_SurfTemp(0:imax-1, 1:jmax) :real(DP), intent(in)
: 地表面温度. Surface temperature
xyz_DelAtmMass(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(in)
: Atmospheric mass of layers
xyz_QH2OVap(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(in)
: Specific humidity
xyr_RadLUwFlux(0:imax-1, 1:jmax, 0:kmax) :real(DP), intent(out)
: 長波フラックス. Upward longwave flux
xyr_RadLDwFlux(0:imax-1, 1:jmax, 0:kmax) :real(DP), intent(out)
: 長波フラックス. Downward longwave flux
xyra_DelRadLUwFlux(0:imax-1, 1:jmax, 0:kmax, 0:1) :real(DP), intent(out)
: 長波地表温度変化.
xyra_DelRadLDwFlux(0:imax-1, 1:jmax, 0:kmax, 0:1) :real(DP), intent(out)
: 長波地表温度変化.

長波フラックスの計算

Calculate long wave flux

[Source]

  subroutine OLD_RadSimpleLWFlux( xy_SurfEmis, xyz_Temp, xyz_Press, xy_SurfTemp, xyz_DelAtmMass, xyz_QH2OVap, xyr_RadLUwFlux, xyr_RadLDwFlux, xyra_DelRadLUwFlux, xyra_DelRadLDwFlux )
    !
    ! 長波フラックスの計算
    !
    ! Calculate long wave flux
    !

    ! モジュール引用 ; USE statements
    !

    ! 物理・数学定数設定
    ! Physical and mathematical constants settings
    !
    use constants0, only: PI, StB                   ! $ \sigma_{SB} $ . 
                              ! ステファンボルツマン定数. 
                              ! Stefan-Boltzmann constant

    ! プランク関数の計算
    ! Calculate Planck function
    !
    use planck_func, only : Integ_PF_GQ_Array3D, Integ_PF_GQ_Array2D, Integ_DPFDT_GQ_Array2D

    ! 散乱を無視した放射伝達方程式
    ! Radiative transfer equation without considering scattering
    !
    use rad_rte_nonscat, only : RadRTENonScat


    ! 宣言文 ; Declaration statements
    !
    real(DP), intent(in):: xy_SurfEmis       (0:imax-1, 1:jmax)
                              ! 地表面射出率.
                              ! Surface emissivity
    real(DP), intent(in):: xyz_Temp          (0:imax-1, 1:jmax, 1:kmax)
                              ! $ T $ .     温度. Temperature
    real(DP), intent(in):: xyz_Press         (0:imax-1, 1:jmax, 1:kmax)
                              ! $ T $ .     温度. Temperature

    real(DP), intent(in):: xy_SurfTemp       (0:imax-1, 1:jmax)
                              ! 地表面温度. 
                              ! Surface temperature
    real(DP), intent(in):: xyz_DelAtmMass    (0:imax-1, 1:jmax, 1:kmax)
                              ! 
                              ! Atmospheric mass of layers
    real(DP), intent(in):: xyz_QH2OVap       (0:imax-1, 1:jmax, 1:kmax)
                              !
                              ! Specific humidity
    real(DP), intent(out):: xyr_RadLUwFlux     (0:imax-1, 1:jmax, 0:kmax)
                              ! 長波フラックス. 
                              ! Upward longwave flux
    real(DP), intent(out):: xyr_RadLDwFlux     (0:imax-1, 1:jmax, 0:kmax)
                              ! 長波フラックス. 
                              ! Downward longwave flux
    real(DP), intent(out):: xyra_DelRadLUwFlux (0:imax-1, 1:jmax, 0:kmax, 0:1)
                              ! 長波地表温度変化. 
                              ! 
    real(DP), intent(out):: xyra_DelRadLDwFlux (0:imax-1, 1:jmax, 0:kmax, 0:1)
                              ! 長波地表温度変化. 
                              ! 

    ! 作業変数
    ! Work variables
    !
    real(DP) :: xyz_DelOptDepDryCom(0:imax-1, 1:jmax, 1:kmax)
    real(DP) :: xyz_DelOptDepH2OVap(0:imax-1, 1:jmax, 1:kmax)
    real(DP) :: xyz_TransEachLayer (0:imax-1, 1:jmax, 1:kmax)
    real(DP) :: xyrr_Trans         (0:imax-1, 1:jmax, 0:kmax, 0:kmax)
                              ! 透過係数. 
                              ! Transmission coefficient
    real(DP) :: xyz_IntPF   (0:imax-1, 1:jmax, 1:kmax)
    real(DP) :: xy_SurfIntPF(0:imax-1, 1:jmax)
    real(DP) :: xy_IntDPFDT0(0:imax-1, 1:jmax)
    real(DP) :: xy_IntDPFDT1(0:imax-1, 1:jmax)

    real(DP) :: xyr_RadUwFlux     (0:imax-1, 1:jmax, 0:kmax)
                              ! 長波フラックス. 
                              ! Upward longwave flux
    real(DP) :: xyr_RadDwFlux     (0:imax-1, 1:jmax, 0:kmax)
                              ! 長波フラックス. 
                              ! Downward longwave flux
    real(DP) :: xyra_DelRadUwFlux (0:imax-1, 1:jmax, 0:kmax, 0:1)
                              ! 長波地表温度変化. 
                              ! 
    real(DP) :: xyra_DelRadDwFlux (0:imax-1, 1:jmax, 0:kmax, 0:1)
                              ! 長波地表温度変化. 
                              ! 

    real(DP) :: WNs
    real(DP) :: WNe

    integer:: k, kk           ! 鉛直方向に回る DO ループ用作業変数
                              ! Work variables for DO loop in vertical direction
    integer:: n               ! 波長について回る DO ループ用作業変数
                              ! Work variables for DO loop in wavenumber bands

    ! 実行文 ; Executable statement
    !


    !   Initialization
    !
    xyr_RadLUwFlux     = 0.0_DP
    xyr_RadLDwFlux     = 0.0_DP
    xyra_DelRadLUwFlux = 0.0_DP
    xyra_DelRadLDwFlux = 0.0_DP
    !
    LOOP_BAND_RTE : do n = 1, nbmax


      ! $ \pi B $, $ \pi DBDT $ の計算
      ! Calculate $ \pi B $ and $ \pi DBDT $
      !
      if ( nbmax == 1 ) then

        xy_SurfIntPF = xy_SurfEmis * StB * ( xy_SurfTemp**4 )
        xyz_IntPF    =               StB * ( xyz_Temp**4 )
        xy_IntDPFDT0 = xy_SurfEmis * 4.0_DP * xy_SurfIntPF / xy_SurfTemp
        xy_IntDPFDT1 =               4.0_DP * xyz_IntPF(:,:,1) / xyz_Temp(:,:,1)

      else

        WNs = a_WNBnds(n-1)
        WNe = a_WNBnds(n  )
        call Integ_PF_GQ_Array3D( WNs, WNe, NumGaussNode, 0, imax-1, 1, jmax, 1, kmax, xyz_Temp, xyz_IntPF )
        call Integ_PF_GQ_Array2D( WNs, WNe, NumGaussNode, 0, imax-1, 1, jmax, xy_SurfTemp, xy_SurfIntPF )
        call Integ_DPFDT_GQ_Array2D( WNs, WNe, NumGaussNode, 0, imax-1, 1, jmax, xyz_Temp(:,:,1), xy_IntDPFDT1 )
        call Integ_DPFDT_GQ_Array2D( WNs, WNe, NumGaussNode, 0, imax-1, 1, jmax, xy_SurfTemp, xy_IntDPFDT0 )

        xy_SurfIntPF = xy_SurfEmis * PI * xy_SurfIntPF
        xyz_IntPF    =               PI * xyz_IntPF
        xy_IntDPFDT0 = xy_SurfEmis * PI * xy_IntDPFDT0
        xy_IntDPFDT1 =               PI * xy_IntDPFDT1

      end if


      ! 光学的厚さの計算
      ! Calculate optical depth
      !
      xyz_DelOptDepDryCom = a_AbsCoefDryCom(n) * ( xyz_Press / a_RefPressDryCom(n) )**a_PressScaleIndDryCom(n) * xyz_DelAtmMass
      xyz_DelOptDepH2OVap = a_AbsCoefH2OVap(n) * ( xyz_Press / a_RefPressH2OVap(n) )**a_PressScaleIndH2OVap(n) * xyz_DelAtmMass * xyz_QH2OVap


      ! 透過関数の計算
      ! Calculate transmission functions
      !
      xyz_TransEachLayer = exp( - DiffFact * ( xyz_DelOptDepDryCom + xyz_DelOptDepH2OVap ) )
      do k = 0, kmax
        do kk = k, k
          xyrr_Trans(:,:,k,kk) = 1.0d0
        end do
        do kk = k+1, kmax
          xyrr_Trans(:,:,k,kk) = xyrr_Trans(:,:,k,kk-1) * xyz_TransEachLayer(:,:,kk)
        end do
        do kk = 0, k-1
          xyrr_Trans(:,:,k,kk) = xyrr_Trans(:,:,kk,k)
        end do
      end do


      call RadRTENonScat( xyz_IntPF, xy_SurfIntPF, xy_IntDPFDT1, xy_IntDPFDT0, xyrr_Trans, xyr_RadUwFlux, xyr_RadDwFlux, xyra_DelRadUwFlux, xyra_DelRadDwFlux )


      xyr_RadLUwFlux     = xyr_RadLUwFlux     + xyr_RadUwFlux
      xyr_RadLDwFlux     = xyr_RadLDwFlux     + xyr_RadDwFlux
      xyra_DelRadLUwFlux = xyra_DelRadLUwFlux + xyra_DelRadUwFlux
      xyra_DelRadLDwFlux = xyra_DelRadLDwFlux + xyra_DelRadDwFlux

    end do LOOP_BAND_RTE


  end subroutine OLD_RadSimpleLWFlux
Old_Flux_saved
Variable :
Old_Flux_saved = .false. :logical, save
: 一度計算したフラックスを保存したことを示すフラグ. Flag for saving of flux calculated once
RadActDryComMMR
Variable :
RadActDryComMMR :real(DP), save
: Mass mixing ratio of radiative active dry component
a_AbsCoefDryCom
Variable :
a_AbsCoefDryCom(1:MaxNmlArySize) :real(DP), save
: $ bar{k}_R $ . 空気の吸収係数 Absorption coefficient of dry component
a_AbsCoefH2OVap
Variable :
a_AbsCoefH2OVap(1:MaxNmlArySize) :real(DP), save
: $ k_R $ . 水の吸収係数 Absorption coefficient of water vapor
a_PressScaleIndDryCom
Variable :
a_PressScaleIndDryCom(1:MaxNmlArySize) :real(DP), save
: Pressure scaling index for dry component
a_PressScaleIndH2OVap
Variable :
a_PressScaleIndH2OVap(1:MaxNmlArySize) :real(DP), save
: Pressure scaling index for water vapor
a_RefPressDryCom
Variable :
a_RefPressDryCom(1:MaxNmlArySize) :real(DP), save
: Reference pressure for dry component.
a_RefPressH2OVap
Variable :
a_RefPressH2OVap(1:MaxNmlArySize) :real(DP), save
: Reference pressure for water vapor
a_WNBnds
Variable :
a_WNBnds(0:MaxNmlArySize) :real(DP), save
: Wavenumber bounds for bands
module_name
Constant :
module_name = ‘rad_simple_LW :character(*), parameter
: モジュールの名称. Module name
nbmax
Variable :
nbmax :integer , save
: 長波バンド数. Number of long wave band
rad_simple_LW_inited
Variable :
rad_simple_LW_inited = .false. :logical, save
: 初期設定フラグ. Initialization flag
version
Constant :
version = ’$Name: $’ // ’$Id: rad_simple_LW.f90,v 1.3 2013/05/25 06:47:33 yot Exp $’ :character(*), parameter
: モジュールのバージョン Module version