Class radiation_band
In: radiation/radiation_band.f90

放射フラックス (バンドモデル)

Radiation flux (band model)

Note that Japanese and English are described in parallel.

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

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

Procedures List

RadiationFlux :放射フラックスの計算
RadiationCorrect :放射フラックス補正
RadiationDTempDt :放射フラックスによる温度変化の計算
———— :————
RadiationFlux :Calculate radiation flux
RadiationCorrect :Radiation flux correction
RadiationDTempDt :Calculate temperature tendency with radiation flux

NAMELIST

NAMELIST#radiation_band_nml

Methods

Included Modules

gridset dc_date_types dc_types namelist_util dc_message constants timeset dc_date dc_trace gtool_historyauto axesset dc_iounit dc_string

Public Instance methods

Subroutine :
xyz_DTempDt(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(in)
: $ DP{T}{t} $ . 温度変化. Temperature tendency
xy_DSurfTempDt(0:imax-1, 1:jmax) :real(DP), intent(in)
: 地表面温度変化率. Surface temperature tendency
xyra_DelRadLFlux(0:imax-1, 1:jmax, 0:kmax, 0:1) :real(DP), intent(in)
: 長波地表温度変化. Surface temperature tendency with longwave
xyr_RadLFlux(0:imax-1, 1:jmax, 0:kmax) :real(DP), intent(inout)
: 長波フラックス. Longwave flux

放射フラックスの補正を行います.

Radiation flux is corrected.

[Source]

  subroutine RadiationCorrect( xyz_DTempDt, xy_DSurfTempDt, xyra_DelRadLFlux, xyr_RadLFlux )
    !
    ! 放射フラックスの補正を行います. 
    !
    ! Radiation flux is corrected. 
    !

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

    ! 時刻管理
    ! Time control
    !
    use timeset, only: DelTime, TimesetClockStart, TimesetClockStop

    ! 宣言文 ; Declaration statements
    !
    implicit none
    real(DP), intent(in):: xyz_DTempDt (0:imax-1, 1:jmax, 1:kmax)
                              ! $ \DP{T}{t} $ . 温度変化. 
                              ! Temperature tendency
    real(DP), intent(in):: xy_DSurfTempDt (0:imax-1, 1:jmax)
                              ! 地表面温度変化率. 
                              ! Surface temperature tendency
    real(DP), intent(in):: xyra_DelRadLFlux (0:imax-1, 1:jmax, 0:kmax,  0:1)
                              ! 長波地表温度変化. 
                              ! Surface temperature tendency with longwave
    real(DP), intent(inout):: xyr_RadLFlux (0:imax-1, 1:jmax, 0:kmax)
                              ! 長波フラックス. 
                              ! Longwave flux

    ! 作業変数
    ! Work variables
    !
    integer:: k               ! 鉛直方向に回る DO ループ用作業変数
                              ! Work variables for DO loop in vertical direction

    ! 実行文 ; Executable statement
    !

    ! 計算時間計測開始
    ! Start measurement of computation time
    !
    call TimesetClockStart( module_name )

    ! 初期化
    ! Initialization
    !
    if ( .not. radiation_band_inited ) call RadiationInit

    ! 放射フラックスの補正
    ! Correct radiation flux
    !
    do k = 0, kmax
      xyr_RadLFlux(:,:,k) = xyr_RadLFlux(:,:,k) + (   xy_DSurfTempDt     * xyra_DelRadLFlux(:,:,k,0) + xyz_DTempDt(:,:,1) * xyra_DelRadLFlux(:,:,k,1) ) * 2.0_DP * DelTime
    end do

    ! 計算時間計測一時停止
    ! Pause measurement of computation time
    !
    call TimesetClockStop( module_name )

  end subroutine RadiationCorrect
Subroutine :
xyr_RadLFlux(0:imax-1, 1:jmax, 0:kmax) :real(DP), intent(in)
: 長波フラックス. Longwave flux
xyr_RadSFlux(0:imax-1, 1:jmax, 0:kmax) :real(DP), intent(in)
: 短波 (日射) フラックス. Shortwave (insolation) flux
xyr_Press(0:imax-1, 1:jmax, 0:kmax) :real(DP), intent(in)
: $ hat{p} $ . 気圧 (半整数レベル). Air pressure (half level)
xyz_DTempDtRadL(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(out)
: 長波加熱率. Temperature tendency with longwave
xyz_DTempDtRadS(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(out)
: 短波加熱率. Temperature tendency with shortwave

放射による温度変化率を計算します.

Temperature tendency with radiation is calculated.

[Source]

  subroutine RadiationDTempDt( xyr_RadLFlux, xyr_RadSFlux, xyr_Press, xyz_DTempDtRadL, xyz_DTempDtRadS )
    !
    ! 放射による温度変化率を計算します. 
    ! 
    ! Temperature tendency with radiation is calculated. 
    !

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

    ! 物理定数設定
    ! Physical constants settings
    !
    use constants, only: Grav, CpDry    ! $ C_p $ [J kg-1 K-1]. 
                                  ! 乾燥大気の定圧比熱. 
                                  ! Specific heat of air at constant pressure

    ! 時刻管理
    ! Time control
    !
    use timeset, only: TimeN, TimesetClockStart, TimesetClockStop

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

    ! 宣言文 ; Declaration statements
    !
    implicit none
    real(DP), intent(in):: xyr_RadLFlux (0:imax-1, 1:jmax, 0:kmax)
                              ! 長波フラックス. 
                              ! Longwave flux
    real(DP), intent(in):: xyr_RadSFlux (0:imax-1, 1:jmax, 0:kmax)
                              ! 短波 (日射) フラックス. 
                              ! Shortwave (insolation) flux
    real(DP), intent(in):: xyr_Press    (0:imax-1, 1:jmax, 0:kmax)
                              ! $ \hat{p} $ . 気圧 (半整数レベル). 
                              ! Air pressure (half level)
    real(DP), intent(out):: xyz_DTempDtRadL (0:imax-1, 1:jmax, 1:kmax)
                              ! 長波加熱率. 
                              ! Temperature tendency with longwave
    real(DP), intent(out):: xyz_DTempDtRadS (0:imax-1, 1:jmax, 1:kmax)
                              ! 短波加熱率. 
                              ! Temperature tendency with shortwave

    ! 作業変数
    ! Work variables
    !
    integer:: k               ! 鉛直方向に回る DO ループ用作業変数
                              ! Work variables for DO loop in vertical direction

    ! 実行文 ; Executable statement
    !

    ! 計算時間計測開始
    ! Start measurement of computation time
    !
    call TimesetClockStart( module_name )

    ! 初期化
    ! Initialization
    !
    if ( .not. radiation_band_inited ) call RadiationInit

    ! 放射冷却率の演算
    ! Calculate radiation cooling rate
    !
    do k = 1, kmax
      xyz_DTempDtRadL(:,:,k) = (     xyr_RadLFlux(:,:,k-1) - xyr_RadLFlux(:,:,k) ) / ( xyr_Press(:,:,k-1)    - xyr_Press(:,:,k) ) / CpDry * Grav

      xyz_DTempDtRadS(:,:,k) = (     xyr_RadSFlux(:,:,k-1) - xyr_RadSFlux(:,:,k) ) / ( xyr_Press(:,:,k-1)    - xyr_Press(:,:,k) ) / CpDry * Grav
    end do

    ! ヒストリデータ出力
    ! History data output
    !
    call HistoryAutoPut( TimeN, 'RadLFlux',    xyr_RadLFlux )
    call HistoryAutoPut( TimeN, 'RadSFlux',    xyr_RadSFlux )
    call HistoryAutoPut( TimeN, 'DTempDtRadL', xyz_DTempDtRadL )
    call HistoryAutoPut( TimeN, 'DTempDtRadS', xyz_DTempDtRadS )

    ! 計算時間計測一時停止
    ! Pause measurement of computation time
    !
    call TimesetClockStop( module_name )

  end subroutine RadiationDTempDt
Subroutine :
xyz_Temp(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(in)
: $ T $ . 温度. Temperature
xyz_QVap(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(in)
: $ q $ . 比湿. Specific humidity
xyr_Press(0:imax-1, 1:jmax, 0:kmax) :real(DP), intent(in)
: $ hat{p} $ . 気圧 (半整数レベル). Air pressure (half level)
xy_SurfTemp(0:imax-1, 1:jmax) :real(DP), intent(in)
: 地表面温度. Surface temperature
xy_SurfAlbedo(0:imax-1, 1:jmax) :real(DP), intent(in)
: 地表アルベド. Surface albedo
xyr_RadLFlux(0:imax-1, 1:jmax, 0:kmax) :real(DP), intent(out)
: 長波フラックス. Longwave flux
xyr_RadSFlux(0:imax-1, 1:jmax, 0:kmax) :real(DP), intent(out)
: 短波 (日射) フラックス. Shortwave (insolation) flux
xya_SurfRadLMtx(0:imax-1, 1:jmax, -1:1) :real(DP), intent(out)
: $ T $ . 陰解行列: 地表. implicit matrix: surface
xyra_DelRadLFlux(0:imax-1, 1:jmax, 0:kmax, 0:1) :real(DP), intent(out)
: 長波地表温度変化. Surface temperature tendency with longwave

温度, 比湿, 気圧から, 放射フラックスを計算します.

Calculate radiation flux from temperature, specific humidity, and air pressure.

[Source]

  subroutine RadiationFlux( xyz_Temp, xyz_QVap, xyr_Press, xy_SurfTemp, xy_SurfAlbedo, xyr_RadLFlux, xyr_RadSFlux, xya_SurfRadLMtx, xyra_DelRadLFlux )
    !
    ! 温度, 比湿, 気圧から, 放射フラックスを計算します. 
    !
    ! Calculate radiation flux from temperature, specific humidity, and 
    ! air pressure. 
    !

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

    ! 物理定数設定
    ! Physical constants settings
    !
    use constants, only: Grav ! $ g $ [m s-2]. 
                              ! 重力加速度. 
                              ! Gravitational acceleration

    ! 時刻管理
    ! Time control
    !
    use timeset, only: TimeN, TimesetClockStart, TimesetClockStop

    ! 日付および時刻の取り扱い
    ! Date and time handler
    !
    use dc_date, only: operator(-), operator(>=), toChar

    ! デバッグ用ユーティリティ
    ! Utilities for debug
    !
    use dc_trace, only: DbgMessage, BeginSub, EndSub, Debug

    ! 宣言文 ; Declaration statements
    !
    implicit none

    real(DP), intent(in):: xyz_Temp  (0:imax-1, 1:jmax, 1:kmax)
                              ! $ T $ .     温度. Temperature
    real(DP), intent(in):: xyz_QVap  (0:imax-1, 1:jmax, 1:kmax)
                              ! $ q $ .     比湿. Specific humidity
    real(DP), intent(in):: xyr_Press (0:imax-1, 1:jmax, 0:kmax)
                              ! $ \hat{p} $ . 気圧 (半整数レベル). 
                              ! Air pressure (half level)
    real(DP), intent(in):: xy_SurfTemp (0:imax-1, 1:jmax)
                              ! 地表面温度. 
                              ! Surface temperature
    real(DP), intent(in):: xy_SurfAlbedo (0:imax-1, 1:jmax)
                              ! 地表アルベド. 
                              ! Surface albedo

    real(DP), intent(out):: xyr_RadLFlux (0:imax-1, 1:jmax, 0:kmax)
                              ! 長波フラックス. 
                              ! Longwave flux
    real(DP), intent(out):: xyr_RadSFlux (0:imax-1, 1:jmax, 0:kmax)
                              ! 短波 (日射) フラックス. 
                              ! Shortwave (insolation) flux
    real(DP), intent(out):: xya_SurfRadLMtx (0:imax-1, 1:jmax, -1:1)
                              ! $ T $ .  陰解行列: 地表. 
                              ! implicit matrix: surface
    real(DP), intent(out):: xyra_DelRadLFlux (0:imax-1, 1:jmax, 0:kmax, 0:1)
                              ! 長波地表温度変化. 
                              ! Surface temperature tendency with longwave

    ! 作業変数
    ! Work variables
    !
    real(DP):: xyr_TauQVap (0:imax-1, 1:jmax, 0:kmax)
                              ! $ \int_z^{\infty} \rho q \, dz $ . 
                              ! 水に関する光学的厚さ 
                              ! $ \tau = \int_z^{\infty} k \rho q \, dz $
                              ! を吸収係数 $ k $ (高度 $ z $ に拠らず一定) 
                              ! で割ったもの.
                              ! 
                              ! Optical depth of water 
                              ! $ \tau = \int_z^{\infty} k \rho q \, dz $
                              ! divided by absorption coefficient $ k $ 
                              ! (this does not depend to height $ z $
                              ! and constant) of water
                              ! 
    real(DP):: xyr_TauDryAir (0:imax-1, 1:jmax, 0:kmax)
                              ! $ \int_z^{\infty} \rho \, dz $ . 
                              ! 空気に関する光学的厚さ 
                              ! $ \tau = \int_z^{\infty} k \rho \, dz $
                              ! を吸収係数 $ k $ (高度 $ z $ に拠らず一定) 
                              ! で割ったもの.
                              ! 
                              ! Optical depth of air 
                              ! $ \tau = \int_z^{\infty} k \rho \, dz $
                              ! divided by absorption coefficient $ k $ 
                              ! (this does not depend to height $ z $
                              ! and constant) of air
                              ! 

    integer:: k               ! 鉛直方向に回る DO ループ用作業変数
                              ! Work variables for DO loop in vertical direction
    character(STRING):: str_debug
                              ! デバッグ用変数
                              ! Variable for debug

    ! 実行文 ; Executable statement
    !

    ! 計算時間計測開始
    ! Start measurement of computation time
    !
    call TimesetClockStart( module_name )

    ! 初期化
    ! Initialization
    !
    if ( .not. radiation_band_inited ) call RadiationInit

    ! 光学的厚さを吸収係数で割ったものの計算
    ! Calculate optical depth divided by absorption coefficient
    !
    xyr_TauQVap   = 0.0_DP
    xyr_TauDryAir = 0.0_DP

    do k = kmax-1, 0, -1
      xyr_TauQVap(:,:,k) = xyr_TauQVap(:,:,k+1) + xyz_QVap(:,:,k+1) * ( xyr_Press(:,:,k) - xyr_Press(:,:,k+1) ) / Grav

      xyr_TauDryAir(:,:,k) = xyr_TauDryAir(:,:,k+1) + ( xyr_Press(:,:,k) - xyr_Press(:,:,k+1) ) / Grav
    end do

    ! 長波フラックスの算出
    ! Calculate long wave flux
    !
    if ( TimeN - PrevTimeLong >= IntTimeLong .or. .not. Old_Flux_saved ) then

      PrevTimeLong = TimeN
      if ( Debug() ) then
        str_debug = toChar(TimeN)
        call DbgMessage( '%c: LongFlux is calculated at %c', c1 = module_name, c2 = trim(str_debug) )
      end if

      call LongFlux( xyz_Temp, xy_SurfTemp, xyr_TauQVap, xyr_TauDryAir, xyr_RadLFlux, xyra_DelRadLFlux )                     ! (out)

    ! 前回の値を利用
    ! Use values in last time
    !
    else

      if ( Debug() ) then
        str_debug = toChar(TimeN)
        call DbgMessage( '%c: LongFlux is not calculated at %c. Save values are used.', c1 = module_name, c2 = trim(str_debug) )
      end if

      xyr_RadLFlux = xyr_RadLFluxSave
      xyra_DelRadLFlux = xyra_DelRadLFluxSave
       
      do k = 0, kmax
        xyr_RadLFlux(:,:,k) = xyr_RadLFlux(:,:,k) + xyra_DelRadLFlux(:,:,k,1) * ( xyz_Temp(:,:,1) - xy_TempSave )

        xyra_DelRadLFlux(:,:,k,1) = xyra_DelRadLFlux(:,:,k,1) / ( xy_TempSave**3 ) * ( xyz_Temp(:,:,1)**3 )
      end do
    end if

    ! 長波陰解用行列の算出
    ! Calculate long wave implicit matrix
    !
    xya_SurfRadLMtx(:,:,0)  = xyra_DelRadLFlux(:,:,0,0)
    xya_SurfRadLMtx(:,:,1)  = xyra_DelRadLFlux(:,:,0,1)
    xya_SurfRadLMtx(:,:,-1) = 0.0_DP


    ! 短波 (日射) フラックスの算出
    ! Calculate short wave (insolation)
    !
    if ( TimeN - PrevTimeShort >= IntTimeShort .or. .not. Old_Flux_saved ) then

      PrevTimeShort = TimeN
      if ( Debug() ) then
        str_debug = toChar(TimeN)
        call DbgMessage( '%c: ShortFlux is calculated at %c', c1 = module_name, c2 = trim(str_debug) )
      end if

      ! 短波フラックスの計算
      ! Calculate short wave (insolation) flux
      !
      xyr_RadSFlux = 0.0_DP
      xyr_RadSFlux(:,:,kmax) = xy_IncomRadSFlux

      call ShortFlux( xyr_TauQVap, xyr_TauDryAir, xy_SurfAlbedo, xyr_RadSFlux )                               ! (inout)

    ! 前回の値を利用
    ! Use values in last time
    ! 
    else

      if ( Debug() ) then
        str_debug = toChar(TimeN)
        call DbgMessage( '%c: ShortFlux is not calculated at %c. Save values are used.', c1 = module_name, c2 = trim(str_debug) )
      end if

      xyr_RadSFlux = xyr_RadSFluxSave
    end if


    ! 今回計算した値を保存
    ! Save calculated values in this time
    !
    xy_TempSave          = xyz_Temp (:,:,1)
    xyr_RadLFluxSave     = xyr_RadLFlux
    xyr_RadSFluxSave     = xyr_RadSFlux
    xyra_DelRadLFluxSave = xyra_DelRadLFlux

    Old_Flux_saved = .true.

    ! 計算時間計測一時停止
    ! Pause measurement of computation time
    !
    call TimesetClockStop( module_name )

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

Private Instance methods

AtmosAlbedo
Variable :
AtmosAlbedo :real(DP), save
: 大気アルベド. Albedo of air
DelTimeLongUnit
Variable :
DelTimeLongUnit :character(STRING), save
: 長波フラックスを計算する時間間隔の単位. Unit of interval step of long wave flux calculation
DelTimeLongValue
Variable :
DelTimeLongValue :real(DP), save
: 長波フラックスを計算する時間間隔の数値. Value of interval step of long wave flux calculation
DelTimeShortUnit
Variable :
DelTimeShortUnit :character(STRING), save
: 短波 (日射) フラックスを計算する時間間隔の単位. Unit of interval step of short wave (insolation) flux calculation
DelTimeShortValue
Variable :
DelTimeShortValue :real(DP), save
: 短波 (日射) フラックスを計算する時間間隔の数値. Value of interval step of short wave (insolation) flux calculation
IncomAIns
Variable :
IncomAIns :real(DP), save
: $ A_{ins} $ . 年平均入射の係数. Coefficient of annual mean incoming radiation.
IncomAZet
Variable :
IncomAZet :real(DP), save
: $ A_{zeta} $ . 年平均入射角の係数. AIns に同じ. Coefficient of annual mean incoming radiation. Same as "AIns".
IncomBIns
Variable :
IncomBIns :real(DP), save
: $ B_{ins} $ . 年平均入射の係数. AIns に同じ. Coefficient of annual mean incoming radiation. Same as "AIns".
IncomBZet
Variable :
IncomBZet :real(DP), save
: $ B_{zeta} $ . 年平均入射角の係数. AIns に同じ. Coefficient of annual mean incoming radiation. Same as "AIns".
Subroutine :

依存モジュールの初期化チェック

Check initialization of dependency modules

[Source]

  subroutine InitCheck
    !
    ! 依存モジュールの初期化チェック
    !
    ! Check initialization of dependency modules

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

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

    ! 格子点設定
    ! Grid points settings
    !
    use gridset, only: gridset_inited

    ! 物理定数設定
    ! Physical constants settings
    !
    use constants, only: constants_inited

    ! 座標データ設定
    ! Axes data settings
    !
    use axesset, only: axesset_inited

    ! 時刻管理
    ! Time control
    !
    use timeset, only: timeset_inited

    ! 実行文 ; Executable statement
    !

    if ( .not. namelist_util_inited ) call MessageNotify( 'E', module_name, '"namelist_util" module is not initialized.' )

    if ( .not. gridset_inited ) call MessageNotify( 'E', module_name, '"gridset" module is not initialized.' )

    if ( .not. constants_inited ) call MessageNotify( 'E', module_name, '"constants" module is not initialized.' )

    if ( .not. axesset_inited ) call MessageNotify( 'E', module_name, '"axesset" module is not initialized.' )

    if ( .not. timeset_inited ) call MessageNotify( 'E', module_name, '"timeset" module is not initialized.' )

  end subroutine InitCheck
IntTimeLong
Variable :
IntTimeLong :type(DC_DIFFTIME), save
: 長波フラックスを計算する時間間隔. Interval time of long wave flux calculation
IntTimeShort
Variable :
IntTimeShort :type(DC_DIFFTIME), save
: 短波フラックスを計算する時間間隔. Interval time of short wave flux calculation
LongAbsorpCoeffDryAir
Variable :
LongAbsorpCoeffDryAir(1:MaxNmlArySize) :real(DP), save
: $ bar{k}_R $ . 空気の吸収係数. Absorption coefficient of air.
LongAbsorpCoeffQVap
Variable :
LongAbsorpCoeffQVap(1:MaxNmlArySize) :real(DP), save
: $ k_R $ . 水の吸収係数. Absorption coefficient of water.
LongBandNum
Variable :
LongBandNum :integer, save
: 長波バンド数. Number of long wave band
LongBandWeight
Variable :
LongBandWeight(1:MaxNmlArySize) :real(DP), save
: バンドウェイト. Band weight.
Subroutine :
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
xyr_TauQVap(0:imax-1, 1:jmax, 0:kmax) :real(DP), intent(in)
: $ int_z^{infty} rho q \, dz $ . 水に関する光学的厚さ $ tau = int_z^{infty} k rho q \, dz $ を吸収係数 $ k $ (高度 $ z $ に拠らず一定) で割ったもの.

Optical depth of water $ tau = int_z^{infty} k rho q \, dz $ divided by absorption coefficient $ k $ (this does not depend to height $ z $ and constant) of water

xyr_TauDryAir(0:imax-1, 1:jmax, 0:kmax) :real(DP), intent(in)
: $ int_z^{infty} rho \, dz $ . 空気に関する光学的厚さ $ tau = int_z^{infty} k rho \, dz $ を吸収係数 $ k $ (高度 $ z $ に拠らず一定) で割ったもの.

Optical depth of air $ tau = int_z^{infty} k rho \, dz $ divided by absorption coefficient $ k $ (this does not depend to height $ z $ and constant) of air

xyr_RadLFlux(0:imax-1, 1:jmax, 0:kmax) :real(DP), intent(out)
: 長波フラックス. Longwave flux
xyra_DelRadLFlux(0:imax-1, 1:jmax, 0:kmax, 0:1) :real(DP), intent(out)
: 長波地表温度変化. Surface temperature tendency with longwave

長波フラックスの計算

Calculate long wave flux

[Source]

  subroutine LongFlux( xyz_Temp, xy_SurfTemp, xyr_TauQVap, xyr_TauDryAir, xyr_RadLFlux, xyra_DelRadLFlux )
    !
    ! 長波フラックスの計算
    !
    ! Calculate long wave flux
    !

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

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

    ! 宣言文 ; Declaration statements
    !
    implicit none
    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):: xyr_TauQVap (0:imax-1, 1:jmax, 0:kmax)
                              ! $ \int_z^{\infty} \rho q \, dz $ . 
                              ! 水に関する光学的厚さ 
                              ! $ \tau = \int_z^{\infty} k \rho q \, dz $
                              ! を吸収係数 $ k $ (高度 $ z $ に拠らず一定) 
                              ! で割ったもの.
                              ! 
                              ! Optical depth of water 
                              ! $ \tau = \int_z^{\infty} k \rho q \, dz $
                              ! divided by absorption coefficient $ k $ 
                              ! (this does not depend to height $ z $
                              ! and constant) of water
                              ! 
    real(DP), intent(in):: xyr_TauDryAir (0:imax-1, 1:jmax, 0:kmax)
                              ! $ \int_z^{\infty} \rho \, dz $ . 
                              ! 空気に関する光学的厚さ 
                              ! $ \tau = \int_z^{\infty} k \rho \, dz $
                              ! を吸収係数 $ k $ (高度 $ z $ に拠らず一定) 
                              ! で割ったもの.
                              ! 
                              ! Optical depth of air 
                              ! $ \tau = \int_z^{\infty} k \rho \, dz $
                              ! divided by absorption coefficient $ k $ 
                              ! (this does not depend to height $ z $
                              ! and constant) of air
                              ! 
    real(DP), intent(out):: xyr_RadLFlux (0:imax-1, 1:jmax, 0:kmax)
                              ! 長波フラックス. 
                              ! Longwave flux
    real(DP), intent(out):: xyra_DelRadLFlux (0:imax-1, 1:jmax, 0:kmax, 0:1)
                              ! 長波地表温度変化. 
                              ! Surface temperature tendency with longwave

    ! 作業変数
    ! Work variables
    !
    real(DP):: xyr_Trans (0:imax-1, 1:jmax, 0:kmax)
                              ! 透過係数. 
                              ! Transmission coefficient
    real(DP):: xyr_Trans1 (0:imax-1, 1:jmax, 0:kmax)
                              ! 1/2 レベルからの透過係数. 
                              ! Transmission coefficient above 1/2 level
    real(DP):: xyr_Trans2 (0:imax-1, 1:jmax, 0:kmax)
                              ! 3/2 レベルからの透過係数. 
                              ! Transmission coefficient above 3/2 level
    real(DP):: xyz_PiB (0:imax-1, 1:jmax, 1:kmax)
                              ! $ \pi B = \sigma T^{4} $
    real(DP):: xy_SurfPiB (0:imax-1, 1:jmax)
                              ! 地表面の $ \pi B $ . 
                              ! $ \pi B $ on surface

    real(DP):: BandWeightSum  ! バンドウェイトの和
                              ! Sum of band weights

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

    ! 実行文 ; Executable statement
    !

    ! バンドウェイトの設定
    ! Configure band weight
    !
    BandWeightSum = 0.0_DP

    do bn = 1, LongBandNum
      BandWeightSum =  BandWeightSum + LongBandWeight(bn)
    end do

    do bn = 1, LongBandNum
      LongBandWeight(bn) = LongBandWeight(bn) / BandWeightSum
    end do

    ! $ \pi B $ の計算
    ! Calculate $ \pi B $
    !
    xyz_PiB    = StB * ( xyz_Temp**4 )
    xy_SurfPiB = StB * ( xy_SurfTemp**4 )

    do k = 0, kmax

      ! 透過関数計算
      ! Calculate transmission functions
      !
      xyr_Trans = 0.0_DP
      
      do bn = 1, LongBandNum
        do kk = 0, kmax
          xyr_Trans(:,:,kk) = xyr_Trans(:,:,kk) + LongBandWeight(bn) * exp( - LongPathLengthFact * (   LongAbsorpCoeffQVap(bn) * abs(   xyr_TauQVap(:,:,kk) - xyr_TauQVap(:,:,k)  ) + LongAbsorpCoeffDryAir(bn) * abs(   xyr_TauDryAir(:,:,kk) - xyr_TauDryAir(:,:,k)  ) ) )
        end do
      end do

      ! 放射フラックス計算
      ! Calculate radiation flux
      !
      xyr_RadLFlux(:,:,k) = xy_SurfPiB * xyr_Trans(:,:,0)

      do kk = 0, kmax-1
        xyr_RadLFlux(:,:,k) = xyr_RadLFlux(:,:,k) - xyz_PiB(:,:,kk+1) * ( xyr_Trans(:,:,kk) - xyr_Trans(:,:,kk+1) )
      end do

      ! 補正項計算用透過関数計算
      ! Calculate transmission functions for correction terms
      !
      xyr_Trans1(:,:,k) = xyr_Trans(:,:,0)
      xyr_Trans2(:,:,k) = xyr_Trans(:,:,1)

    end do

    ! 長波地表温度変化の計算
    ! Calclate surface temperature tendency with long wave
    !
    do k = 0, kmax
      xyra_DelRadLFlux(:,:,k,0) = 4.0_DP * xy_SurfPiB / xy_SurfTemp * xyr_Trans1(:,:,k)

      xyra_DelRadLFlux(:,:,k,1) = 4.0_DP * xyz_PiB(:,:,1) / xyz_Temp(:,:,1) * ( xyr_Trans2(:,:,k) - xyr_Trans1(:,:,k) )
    end do

  end subroutine LongFlux
LongPathLengthFact
Variable :
LongPathLengthFact :real(DP), save
: 光路長のファクタ. Factor of optical length
Old_Flux_saved
Variable :
Old_Flux_saved = .false. :logical, save
: 一度計算したフラックスを保存したことを示すフラグ. Flug for saving of flux calculated once
PrevTimeLong
Variable :
PrevTimeLong :type(DC_DIFFTIME), save
: 前回長波フラックスを計算した時刻 Time when long wave flux is calculated
PrevTimeShort
Variable :
PrevTimeShort :type(DC_DIFFTIME), save
: 前回短波フラックスを計算した時刻 Time when short wave flux is calculated
Subroutine :

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

"radiation_band" module is initialized. "NAMELIST#radiation_band_nml" is loaded in this procedure.

This procedure input/output NAMELIST#radiation_band_nml .

[Source]

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

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

    ! 時刻管理
    ! Time control
    !
    use timeset, only: StartTime             ! 計算開始時刻. 

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

    ! 日付および時刻の取り扱い
    ! Date and time handler
    !
    use dc_date, only: DCDiffTimeCreate

    ! ファイル入出力補助
    ! 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: StoA

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

    ! 宣言文 ; Declaration statements
    !
    implicit none

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

    ! NAMELIST 変数群
    ! NAMELIST group name
    !
    namelist /radiation_band_nml/ DelTimeLongValue, DelTimeLongUnit, DelTimeShortValue, DelTimeShortUnit, LongBandNum, LongAbsorpCoeffQVap, LongAbsorpCoeffDryAir, LongBandWeight, LongPathLengthFact, ShortBandNum, ShortAbsorpCoeffQVap, ShortAbsorpCoeffDryAir, ShortBandWeight, ShortSecScat, SolarCoeff, AtmosAlbedo, IncomAIns, IncomBIns, IncomAZet, IncomBZet
          !
          ! デフォルト値については初期化手続 "radiation_band#RadiationInit" 
          ! のソースコードを参照のこと. 
          !
          ! Refer to source codes in the initialization procedure
          ! "radiation_band#RadiationInit" for the default values. 
          !

    ! 実行文 ; Executable statement
    !

    if ( radiation_band_inited ) return
    call InitCheck

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

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

    DelTimeLongValue = 3.0_DP
    DelTimeLongUnit  = 'hrs.'

    LongBandNum      = 4
    LongAbsorpCoeffQVap    = -999.9_DP
    LongAbsorpCoeffDryAir  = -999.9_DP
    LongBandWeight         = -999.9_DP
    LongAbsorpCoeffQVap   (1:LongBandNum) = (/ 8.0_DP, 1.0_DP, 0.1_DP, 0.0_DP /)
    LongAbsorpCoeffDryAir (1:LongBandNum) = (/ 0.0_DP, 0.0_DP, 0.0_DP, 5.0e-5_DP /)
    LongBandWeight        (1:LongBandNum) = (/ 0.2_DP, 0.1_DP, 0.1_DP, 0.6_DP /)
    LongPathLengthFact = 1.5_DP

    ! 短波フラックス用情報
    ! Information for short wave flux
    !
    PrevTimeShort = StartTime

    DelTimeShortValue = 1.0_DP
    DelTimeShortUnit  = 'hrs.'

    ShortBandNum = 1
    ShortAbsorpCoeffQVap   = -999.9_DP
    ShortAbsorpCoeffDryAir = -999.9_DP
    ShortBandWeight        = -999.9_DP
    ShortAbsorpCoeffQVap   (1:ShortBandNum) = (/ 0.002_DP /)
    ShortAbsorpCoeffDryAir (1:ShortBandNum) = (/ 0.0_DP /)
    ShortBandWeight        (1:ShortBandNum) = (/ 1.0_DP /)
    ShortSecScat = 1.66_DP

    ! 短波入射用情報
    ! Information for short wave incomming
    !
    SolarCoeff  = 1380.0_DP
    AtmosAlbedo = 0.2_DP
!!$    EpsOrb = 23.5_DP
!!$    EqnOrb = 110.0_DP
    IncomAIns   = 0.127_DP
    IncomBIns   = 0.183_DP
    IncomAZet   = 0.410_DP
    IncomBZet   = 0.590_DP


    ! 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 = radiation_band_nml, iostat = iostat_nml )        ! (out)
      close( unit_nml )

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

    ! 時間間隔の処理
    ! Handle interval time
    !
    call DCDiffTimeCreate( IntTimeLong, DelTimeLongValue, DelTimeLongUnit )   ! (in)
    call DCDiffTimeCreate( IntTimeShort, DelTimeShortValue, DelTimeShortUnit ) ! (in)

    ! バンド数, 吸収係数, バンドウェイトのチェック
    ! Check number of band, absorption coefficients, band weight
    !
    call NmlutilAryValid( module_name, LongAbsorpCoeffQVap, 'LongAbsorpCoeffQVap', LongBandNum,         'LongBandNum' )          ! (in)

    call NmlutilAryValid( module_name, LongAbsorpCoeffDryAir, 'LongAbsorpCoeffDryAir', LongBandNum,           'LongBandNum' )            ! (in)

    call NmlutilAryValid( module_name, LongBandWeight, 'LongBandWeight', LongBandNum,    'LongBandNum' )      ! (in)

    call NmlutilAryValid( module_name, ShortAbsorpCoeffQVap, 'ShortAbsorpCoeffQVap', ShortBandNum,         'ShortBandNum' )          ! (in)

    call NmlutilAryValid( module_name, ShortAbsorpCoeffDryAir, 'ShortAbsorpCoeffDryAir', ShortBandNum,           'ShortBandNum' )            ! (in)

    call NmlutilAryValid( module_name, ShortBandWeight, 'ShortBandWeight', ShortBandNum,    'ShortBandNum' )      ! (in)

    ! 短波入射の計算
    ! Calculate short wave (insolation) incoming radiation
    !
    allocate( xy_IncomRadSFlux (0:imax-1, 1:jmax) )
    allocate( xy_InAngle (0:imax-1, 1:jmax) )

    call ShortIncoming( xy_IncomRadSFlux, xy_InAngle ) ! (out)

    ! 保存用の変数の割り付け
    ! Allocate variables for saving
    !
    allocate( xy_TempSave          (0:imax-1, 1:jmax) )
    allocate( xyr_RadLFluxSave     (0:imax-1, 1:jmax, 0:kmax) )
    allocate( xyr_RadSFluxSave     (0:imax-1, 1:jmax, 0:kmax) )
    allocate( xyra_DelRadLFluxSave (0:imax-1, 1:jmax, 0:kmax, 0:1) )

    ! ヒストリデータ出力のためのへの変数登録
    ! Register of variables for history data output
    !
    call HistoryAutoAddVariable( 'RadLFlux' , (/ 'lon ', 'lat ', 'sigm', 'time'/), 'longwave flux', 'W m-2' )
    call HistoryAutoAddVariable( 'RadSFlux' , (/ 'lon ', 'lat ', 'sigm', 'time' /), 'shortwave flux', 'W m-2' )
    call HistoryAutoAddVariable( 'DTempDtRadL' , (/ 'lon ', 'lat ', 'sig ', 'time' /), 'longwave heating', 'K s-1' )
    call HistoryAutoAddVariable( 'DTempDtRadS' , (/ 'lon ', 'lat ', 'sig ', 'time' /), 'shortwave heating', 'K s-1' )

    ! 印字 ; Print
    !
    call MessageNotify( 'M', module_name, '----- Initialization Messages -----' )
    call MessageNotify( 'M', module_name, 'DelTime:' )
    call MessageNotify( 'M', module_name, '  DelTimeLong  = %f [%c]', d = (/ DelTimeLongValue /), c1 = trim( DelTimeLongUnit ) )
    call MessageNotify( 'M', module_name, '  DelTimeShort = %f [%c]', d = (/ DelTimeShortValue /), c1 = trim( DelTimeShortUnit ) )
!
    call MessageNotify( 'M', module_name, 'LongFlux:' )
    call MessageNotify( 'M', module_name, '  LongBandNum            = %d', i = (/ LongBandNum /) )
    call MessageNotify( 'M', module_name, '  LongAbsorpCoeffQVap    = (/ %*r /)', r = real( LongAbsorpCoeffQVap(1:LongBandNum) ), n = (/ LongBandNum /) )
    call MessageNotify( 'M', module_name, '  LongAbsorpCoeffDryAir  = (/ %*r /)', r = real( LongAbsorpCoeffDryAir(1:LongBandNum) ), n = (/ LongBandNum /) )
    call MessageNotify( 'M', module_name, '  LongBandWeight         = (/ %*r /)', r = real( LongBandWeight(1:LongBandNum) ), n = (/ LongBandNum /) )
    call MessageNotify( 'M', module_name, '  LongPathLengthFact     = %f', d = (/ LongPathLengthFact /) )
!
    call MessageNotify( 'M', module_name, 'ShortFlux:' )
    call MessageNotify( 'M', module_name, '  ShortBandNum           = %d', i = (/ ShortBandNum /) )
    call MessageNotify( 'M', module_name, '  ShortAbsorpCoeffQVap   = (/ %*r /)', r = real( ShortAbsorpCoeffQVap(1:ShortBandNum) ), n = (/ ShortBandNum /) )
    call MessageNotify( 'M', module_name, '  ShortAbsorpCoeffDryAir = (/ %*r /)', r = real( ShortAbsorpCoeffDryAir(1:ShortBandNum) ), n = (/ ShortBandNum /) )
    call MessageNotify( 'M', module_name, '  ShortBandWeight        = (/ %*r /)', r = real( ShortBandWeight(1:ShortBandNum) ), n = (/ ShortBandNum /) )
    call MessageNotify( 'M', module_name, '  ShortSecScat           = %f', d = (/ ShortSecScat /) )
!
    call MessageNotify( 'M', module_name, 'ShortIncomming:' )
    call MessageNotify( 'M', module_name, '  SolarCoeff  = %f', d = (/ SolarCoeff  /) )
    call MessageNotify( 'M', module_name, '  AtmosAlbedo = %f', d = (/ AtmosAlbedo /) )
!!$    call MessageNotify( 'M', module_name, '  EpsOrb      = %f', d = (/ EpsOrb      /) )
!!$    call MessageNotify( 'M', module_name, '  EqnOrb      = %f', d = (/ EqnOrb      /) )
    call MessageNotify( 'M', module_name, '  IncomAIns   = %f', d = (/ IncomAIns   /) )
    call MessageNotify( 'M', module_name, '  IncomBIns   = %f', d = (/ IncomBIns   /) )
    call MessageNotify( 'M', module_name, '  IncomAZet   = %f', d = (/ IncomAZet   /) )
    call MessageNotify( 'M', module_name, '  IncomBZet   = %f', d = (/ IncomBZet   /) )

    call MessageNotify( 'M', module_name, '-- version = %c', c1 = trim(version) )

    radiation_band_inited = .true.
  end subroutine RadiationInit
ShortAbsorpCoeffDryAir
Variable :
ShortAbsorpCoeffDryAir(1:MaxNmlArySize) :real(DP), save
: $ bar{k}_S $ . 空気の吸収係数. Absorption coefficient of air.
ShortAbsorpCoeffQVap
Variable :
ShortAbsorpCoeffQVap(1:MaxNmlArySize) :real(DP), save
: $ k_S $ . 水の吸収係数. Absorption coefficient of water.
ShortBandNum
Variable :
ShortBandNum :integer, save
: 短波バンド数. Number of short wave band
ShortBandWeight
Variable :
ShortBandWeight(1:MaxNmlArySize) :real(DP), save
: バンドウェイト. Band weight.
Subroutine :
xyr_TauQVap(0:imax-1, 1:jmax, 0:kmax) :real(DP), intent(in)
: $ int_z^{infty} rho q \, dz $ . 水に関する光学的厚さ $ tau = int_z^{infty} k rho q \, dz $ を吸収係数 $ k $ (高度 $ z $ に拠らず一定) で割ったもの.

Optical depth of water $ tau = int_z^{infty} k rho q \, dz $ divided by absorption coefficient $ k $ (this does not depend to height $ z $ and constant) of water

xyr_TauDryAir(0:imax-1, 1:jmax, 0:kmax) :real(DP), intent(in)
: $ int_z^{infty} rho \, dz $ . 空気に関する光学的厚さ $ tau = int_z^{infty} k rho \, dz $ を吸収係数 $ k $ (高度 $ z $ に拠らず一定) で割ったもの.

Optical depth of air $ tau = int_z^{infty} k rho \, dz $ divided by absorption coefficient $ k $ (this does not depend to height $ z $ and constant) of air

xy_SurfAlbedo(0:imax-1, 1:jmax) :real(DP), intent(in)
: 地表アルベド. Surface albedo
xyr_RadSFlux(0:imax-1, 1:jmax, 0:kmax) :real(DP), intent(inout)
: 短波 (日射) フラックス. Shortwave (insolation) flux

短波フラックスを計算します.

Calculate short wave flux.

[Source]

  subroutine ShortFlux( xyr_TauQVap, xyr_TauDryAir, xy_SurfAlbedo, xyr_RadSFlux )
    !
    ! 短波フラックスを計算します.
    !
    ! Calculate short wave flux. 
    !

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

    ! 宣言文 ; Declaration statements
    !
    implicit none

    real(DP), intent(in):: xyr_TauQVap (0:imax-1, 1:jmax, 0:kmax)
                              ! $ \int_z^{\infty} \rho q \, dz $ . 
                              ! 水に関する光学的厚さ 
                              ! $ \tau = \int_z^{\infty} k \rho q \, dz $
                              ! を吸収係数 $ k $ (高度 $ z $ に拠らず一定) 
                              ! で割ったもの.
                              ! 
                              ! Optical depth of water 
                              ! $ \tau = \int_z^{\infty} k \rho q \, dz $
                              ! divided by absorption coefficient $ k $ 
                              ! (this does not depend to height $ z $
                              ! and constant) of water
                              ! 
    real(DP), intent(in):: xyr_TauDryAir (0:imax-1, 1:jmax, 0:kmax)
                              ! $ \int_z^{\infty} \rho \, dz $ . 
                              ! 空気に関する光学的厚さ 
                              ! $ \tau = \int_z^{\infty} k \rho \, dz $
                              ! を吸収係数 $ k $ (高度 $ z $ に拠らず一定) 
                              ! で割ったもの.
                              ! 
                              ! Optical depth of air 
                              ! $ \tau = \int_z^{\infty} k \rho \, dz $
                              ! divided by absorption coefficient $ k $ 
                              ! (this does not depend to height $ z $
                              ! and constant) of air
                              ! 
    real(DP), intent(in):: xy_SurfAlbedo (0:imax-1, 1:jmax)
                              ! 地表アルベド. 
                              ! Surface albedo
    real(DP), intent(inout):: xyr_RadSFlux (0:imax-1, 1:jmax, 0:kmax)
                              ! 短波 (日射) フラックス. 
                              ! Shortwave (insolation) flux

    ! 作業変数
    ! Work variables
    !
    real(DP):: BandWeightSum  ! バンドウェイトの和
                              ! Sum of band weights
    integer:: k               ! 鉛直方向に回る DO ループ用作業変数
                              ! Work variables for DO loop in vertical direction
    integer:: bn              ! 波長について回る DO ループ用作業変数
                              ! Work variables for DO loop in wavenumber bands

    ! 実行文 ; Executable statement
    !

    ! バンドウェイトの設定
    ! Configure band weight
    !
    BandWeightSum = 0.0_DP

    do bn = 1, ShortBandNum
      BandWeightSum = BandWeightSum + ShortBandWeight(bn)
    end do

    do bn = 1, ShortBandNum
      ShortBandWeight(bn) = ShortBandWeight(bn) / BandWeightSum
    end do

    do bn = 1, ShortBandNum
      do k = 0, kmax

        ! 各レベルでの下向き透過
        ! Downward transmission on each level
        !
        if ( k /= kmax ) then
          xyr_RadSFlux(:,:,k) = xyr_RadSFlux(:,:,k) + ShortBandWeight(bn) * xyr_RadSFlux(:,:,kmax) * exp( - xy_InAngle(:,:) * (   ShortAbsorpCoeffQVap(bn) * xyr_TauQVap(:,:,k) + ShortAbsorpCoeffDryAir(bn) * xyr_TauDryAir(:,:,k) ) )
        end if
        
        ! 各レベルでの上向き透過
        ! Upward transmission on each level
        !
        xyr_RadSFlux(:,:,k) = xyr_RadSFlux(:,:,k) - ShortBandWeight(bn) * xyr_RadSFlux(:,:,kmax) * exp( - xy_InAngle(:,:) * (   ShortAbsorpCoeffQVap(bn) * xyr_TauQVap(:,:,0) + ShortAbsorpCoeffDryAir(bn) * xyr_TauDryAir(:,:,0) ) ) * xy_SurfAlbedo * exp( - ShortSecScat * (   ShortAbsorpCoeffQVap(bn) * ( xyr_TauQVap(:,:,0) - xyr_TauQVap(:,:,k) ) + ShortAbsorpCoeffDryAir(bn) * ( xyr_TauDryAir(:,:,0) - xyr_TauDryAir(:,:,k) ) ) )
      end do
    end do

    ! 吸収なしの場合
    ! In the case of no absorption
    !
    if ( ShortBandNum == 0 ) then
      do k = 0, kmax
        xyr_RadSFlux(:,:,k) = ( 1.0_DP - xy_SurfAlbedo ) * xyr_RadSFlux(:,:,kmax)
      end do
    end if

  end subroutine ShortFlux
Subroutine :
xy_IncomRadSFlux(0:imax-1, 1:jmax) :real(DP), intent(out)
: 短波 (日射) フラックス. Short wave (insolation) flux
xy_InAngle(0:imax-1, 1:jmax) :real(DP), intent(out)
: sec (入射角). sec (angle of incidence)

短波入射を計算します.

Calculate short wave (insolation) incoming radiation.

[Source]

  subroutine ShortIncoming( xy_IncomRadSFlux, xy_InAngle )
    !
    ! 短波入射を計算します.
    !
    ! Calculate short wave (insolation) incoming radiation. 
    !

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

    ! 座標データ設定
    ! Axes data settings
    !
    use axesset, only: y_Lat  ! $ \varphi $ [rad.] . 緯度. Latitude

    ! 宣言文 ; Declaration statements
    !
    implicit none
    real(DP), intent(out):: xy_IncomRadSFlux (0:imax-1, 1:jmax)
                              ! 短波 (日射) フラックス. 
                              ! Short wave (insolation) flux
    real(DP), intent(out):: xy_InAngle (0:imax-1, 1:jmax)
                              ! sec (入射角). 
                              ! sec (angle of incidence)

    ! 作業変数
    ! Work variables
    !
    integer:: i               ! 経度方向に回る DO ループ用作業変数
                              ! Work variables for DO loop in longitude
    integer:: j               ! 緯度方向に回る DO ループ用作業変数
                              ! Work variables for DO loop in latitude

    ! 実行文 ; Executable statement
    !

    ! 年, 日平均日射の計算
    ! Calculate annual mean, daily mean insolation
    !
    do i = 0, imax - 1
      do j = 1, jmax
        xy_IncomRadSFlux(i,j) = - SolarCoeff * (1.0_DP - AtmosAlbedo ) * ( IncomAIns + IncomBIns * cos( y_Lat(j) )**2 )

        if ( xy_IncomRadSFlux(i,j) < 0.0_DP ) then
          xy_InAngle(i,j) = 1.0_DP / ( IncomAZet + IncomBZet * cos( y_Lat(j) )**2 )
        else
          xy_IncomRadSFlux(i,j) = 0.0_DP
          xy_InAngle(i,j) = 0.0_DP
        end if
      end do
    end do


    ! 季節変化, 日変化がある場合の計算
    ! Calculate with seasonal change and diurnal change
    !

    ! AGCM5 から移植予定. 
    ! Importation from AGCM5 is planed

  end subroutine ShortIncoming
ShortSecScat
Variable :
ShortSecScat :real(DP), save
: 散乱の $ sec zeta $ . $ sec zeta $ of scattering
SolarCoeff
Variable :
SolarCoeff :real(DP), save
: 太陽定数. Solar constant
module_name
Constant :
module_name = ‘radiation_band :character(*), parameter
: モジュールの名称. Module name
version
Constant :
version = ’$Name: dcpam5-20081109-1 $’ // ’$Id: radiation_band.f90,v 1.10 2008-10-07 09:20:27 morikawa Exp $’ :character(*), parameter
: モジュールのバージョン Module version
xy_InAngle
Variable :
xy_InAngle(:,:) :real(DP), allocatable, save
: sec (入射角). sec (angle of incidence)
xy_IncomRadSFlux
Variable :
xy_IncomRadSFlux(:,:) :real(DP), allocatable, save
: 短波 (日射) フラックス. Short wave (insolation) flux
xy_TempSave
Variable :
xy_TempSave(:,:) :real(DP), allocatable, save
: $ T $ . 温度 (保存用). Temperature (for save)
xyr_RadLFluxSave
Variable :
xyr_RadLFluxSave(:,:,:) :real(DP), allocatable, save
: 長波フラックス (保存用). Long wave flux (for save)
xyr_RadSFluxSave
Variable :
xyr_RadSFluxSave(:,:,:) :real(DP), allocatable, save
: 短波 (日射) フラックス (保存用). Short wave (insolation) flux (for save)
xyra_DelRadLFluxSave
Variable :
xyra_DelRadLFluxSave(:,:,:,:) :real(DP), allocatable, save
: 長波地表温度変化 (保存用). Surface temperature tendency with long wave (for save)

[Validate]