Class radiation_short_income
In: radiation/radiation_short_income.f90

短波入射 (太陽入射)

Short wave (insolation) incoming

Note that Japanese and English are described in parallel.

短波入射 (太陽入射) を計算します.

Calculate short wave (insolation) incoming radiation.

Procedures List

ShortIncoming :短波入射 (太陽入射) の計算
———— :————
ShortIncoming :Calculate short wave (insolation) incoming radiation.

NAMELIST

NAMELIST#radiation_short_income_nml

Methods

Included Modules

gridset dc_date_types dc_types namelist_util dc_message axesset dc_iounit gtool_historyauto constants timeset

Public Instance methods

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
    !

    ! 初期化
    ! Initialization
    !
    if ( .not. radiation_short_income_inited ) call ShtIncomeInit

    ! 年, 日平均日射の計算
    ! Calculate annual mean, daily mean insolation
    !
    do i = 0, imax - 1
      do j = 1, jmax
        xy_IncomRadSFlux(i,j) = - SolarConst * (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
radiation_short_income_inited
Variable :
radiation_short_income_inited = .false. :logical, save, public
: 初期設定フラグ. Initialization flag

Private Instance methods

AtmosAlbedo
Variable :
AtmosAlbedo :real(DP), save
: 大気アルベド. Albedo of air
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
Subroutine :

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

"radiation_short_income" module is initialized. "NAMELIST#radiation_short_income_nml" is loaded in this procedure.

This procedure input/output NAMELIST#radiation_short_income_nml .

[Source]

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

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

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

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

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

    ! ヒストリデータ出力
    ! 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_short_income_nml/ SolarConst, AtmosAlbedo, IncomAIns, IncomBIns, IncomAZet, IncomBZet
          !
          ! デフォルト値については初期化手続 "radiation_short_income#ShtIncomeInit" 
          ! のソースコードを参照のこと. 
          !
          ! Refer to source codes in the initialization procedure
          ! "radiation_short_income#ShtIncomeInit" for the default values. 
          !

    ! 実行文 ; Executable statement
    !

    if ( radiation_short_income_inited ) return
    call InitCheck

    ! デフォルト値の設定
    ! Default values settings
    !
    SolarConst  = 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_short_income_nml, iostat = iostat_nml )        ! (out)
      close( unit_nml )

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

    ! 保存用の変数の割り付け
    ! Allocate variables for saving
    !

    ! ヒストリデータ出力のためのへの変数登録
    ! Register of variables for history data output
    !
!!$    call HistoryAutoAddVariable( 'xxxxx' , &
!!$      & (/ 'lon ', 'lat ', 'sig', 'time'/), &
!!$      & 'xxxx', 'W m-2' )

    ! 印字 ; Print
    !
    call MessageNotify( 'M', module_name, '----- Initialization Messages -----' )
    call MessageNotify( 'M', module_name, 'ShortIncomming:' )
    call MessageNotify( 'M', module_name, '  SolarConst  = %f', d = (/ SolarConst  /) )
    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_short_income_inited = .true.
  end subroutine ShtIncomeInit
SolarConst
Variable :
SolarConst :real(DP), save
: 太陽定数. Solar constant
module_name
Constant :
module_name = ‘radiation_short_income :character(*), parameter
: モジュールの名称. Module name
version
Constant :
version = ’$Name: dcpam5-20090218-1 $’ // ’$Id: radiation_short_income.f90,v 1.2 2008-11-16 14:32:22 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

[Validate]