| Class | rad_simple_SW_V2_0 | 
| In: | 
                
                radiation/rad_simple_SW_V2_0.f90
                
         | 
| Subroutine : | |||
| xy_SurfAlbedo(0:imax-1, 1:jmax ) : | real(DP), intent(in ) | ||
| xyz_Press(0:imax-1, 1:jmax, 1:kmax) : | real(DP), intent(in )
  | ||
| xyz_DelAtmMass(0:imax-1, 1:jmax, 1:kmax) : | real(DP), intent(in ) | ||
| xyz_QH2OVap(0:imax-1, 1:jmax, 1:kmax) : | real(DP), intent(in )
  | ||
| xyr_RadSUwFlux(0:imax-1, 1:jmax, 0:kmax) : | real(DP), intent(out) | ||
| xyr_RadSDwFlux(0:imax-1, 1:jmax, 0:kmax) : | real(DP), intent(out) | 
  subroutine RadSimpleSWV20Flux( xy_SurfAlbedo, xyz_Press, xyz_DelAtmMass, xyz_QH2OVap, xyr_RadSUwFlux, xyr_RadSDwFlux )
    ! USE statements
    !
    ! 太陽放射フラックスの設定
    ! Set solar constant
    !
    use set_solarconst, only : SetSolarConst
    ! 短波入射 (太陽入射)
    ! Short wave (insolation) incoming
    !
    use rad_short_income, only : RadShortIncome
    !
    ! Solve radiative transfer equation in two stream approximation
    !
    use rad_rte_two_stream_app, only: RadRTETwoStreamAppSW, RadRTETwoStreamAppHomogAtm
    ! Chou and Lee (1996) による短波放射モデル
    ! Short wave radiation model described by Chou and Lee (1996)
    !
    use rad_CL1996, only : RadCL1996NumBands      , RadCL1996UVVISParams
    real(DP), intent(in ) :: xy_SurfAlbedo    (0:imax-1, 1:jmax )
    real(DP), intent(in ) :: xyz_Press        (0:imax-1, 1:jmax, 1:kmax)
                                ! $ P $ .     圧力. Pressure
    real(DP), intent(in ) :: xyz_DelAtmMass   (0:imax-1, 1:jmax, 1:kmax)
    real(DP), intent(in ) :: xyz_QH2OVap      (0:imax-1, 1:jmax, 1:kmax)
                                !
                                ! Specific humidity
    real(DP), intent(out) :: xyr_RadSUwFlux   (0:imax-1, 1:jmax, 0:kmax)
    real(DP), intent(out) :: xyr_RadSDwFlux   (0:imax-1, 1:jmax, 0:kmax)
    real(DP) :: SolarConst
    real(DP) :: xy_InAngle    (0:imax-1, 1:jmax)
                              ! sec (入射角).
                              ! sec (angle of incidence)
    real(DP) :: DistFromStarScld
                               ! Distance between the central star and the planet
    real(DP) :: DiurnalMeanFactor
    integer  :: nbands1
    integer  :: nbands2
    real(DP) :: UVVISFracSolarFlux
    real(DP) :: UVVISO3AbsCoef
    real(DP) :: UVVISRayScatCoef
    real(DP) :: SolarFluxTOA
!!$    real(DP) :: xyz_SSA       (0:imax-1, 1:jmax, 1:kmax)
!!$    real(DP) :: xyz_AF        (0:imax-1, 1:jmax, 1:kmax)
    real(DP), parameter :: RayScatSinScatAlb = 1.0d0 - 1.0d-10
    real(DP), parameter :: RayScatAsymFact   = 0.0d0
    real(DP)            :: xyz_RayScatDelOptDep( 0:imax-1, 1:jmax, 1:kmax )
    real(DP)            :: xyz_DelOptDepH2OVap(0:imax-1, 1:jmax, 1:kmax)
    real(DP)            :: xyz_DelTotOptDep( 0:imax-1, 1:jmax, 1:kmax )
    real(DP)            :: xyr_TotOptDep   ( 0:imax-1, 1:jmax, 0:kmax )
    real(DP)            :: xyr_RadFlux     ( 0:imax-1, 1:jmax, 0:kmax )
    real(DP)            :: xyr_RadUwFlux   ( 0:imax-1, 1:jmax, 0:kmax )
    real(DP)            :: xyr_RadDwFlux   ( 0:imax-1, 1:jmax, 0:kmax )
    real(DP) :: TotUVVISFracSolarFlux
    integer  :: i
    integer  :: j
    integer  :: k
    integer  :: l
    integer  :: n
    ! 初期化確認
    ! Initialization check
    !
    if ( .not. rad_simple_SW_V2_0_inited ) then
      call MessageNotify( 'E', module_name, 'This module has not been initialized.' )
    end if
    ! 計算時間計測開始
    ! Start measurement of computation time
    !
    call TimesetClockStart( module_name )
    ! 太陽放射フラックスの設定
    ! Set solar constant
    !
    call SetSolarConst( SolarConst )
    ! 短波入射の計算
    ! Calculate short wave (insolation) incoming radiation
    !
    call RadShortIncome( xy_InAngle        = xy_InAngle, DistFromStarScld  = DistFromStarScld, DiurnalMeanFactor = DiurnalMeanFactor )
    call RadCL1996NumBands( nbands1, nbands2 )
    ! Initialization
    !
    xyr_RadSUwFlux = 0.0_DP
    xyr_RadSDwFlux = 0.0_DP
    ! * 14286 to 57143 cm-1 (0.175 to 0.70 micron): 
    !   * Rayleigh scattering, 
    !
    do l = 1, nbands1
      ! UV and Visible optical properties and solar flux
      !
      call RadCL1996UVVISParams( l, UVVISFracSolarFlux, UVVISO3AbsCoef, UVVISRayScatCoef )
      ! Rayleigh scattering coefficient is multiplied by a ratio
      UVVISRayScatCoef = UVVISRayScatCoef * RayleighScatCoefRatio
      SolarFluxTOA = UVVISFracSolarFlux * SolarConst / DistFromStarScld**2 * DiurnalMeanFactor
      ! Rayleigh scattering
      !
      if ( FlagRayleighScattering ) then
        xyz_RayScatDelOptDep = UVVISRayScatCoef * xyz_DelAtmMass
      else
        xyz_RayScatDelOptDep = 0.0d0
      end if
      ! Total optical parameter
      !
      xyz_DelTotOptDep = + xyz_RayScatDelOptDep
      !
      xyr_TotOptDep(:,:,kmax) = 0.0d0
      do k = kmax-1, 0, -1
        xyr_TotOptDep(:,:,k) = xyr_TotOptDep(:,:,k+1) + xyz_DelTotOptDep(:,:,k+1)
      end do
      !
!!$      xyz_SSA = RayScatSinScatAlb
!!$      do k = 1, kmax
!!$        do j = 1, jmax
!!$          do i = 0, imax-1
!!$            if ( xyz_SSA(i,j,k) >= 1.0d0 ) then
!!$              xyz_SSA(i,j,k) = 1.0d0 - 1.0d-10
!!$            end if
!!$          end do
!!$        end do
!!$      end do
!!$      xyz_AF  = RayScatAsymFact
!!$
!!$      call RadRTETwoStreamAppSW(        &
!!$        & xyz_SSA, xyz_AF,              & ! (in)
!!$        & xyr_TotOptDep,                & ! (in)
!!$        & xy_SurfAlbedo,                & ! (in)
!!$        & SolarFluxTOA, xy_InAngle,     & ! (in)
!!$        & xyr_RadUwFlux, xyr_RadDwFlux  & ! (out)
!!$        & )
      call RadRTETwoStreamAppHomogAtm( xy_SurfAlbedo, SolarFluxTOA, xy_InAngle, RayScatSinScatAlb, RayScatAsymFact, xyr_TotOptDep, xyr_RadUwFlux, xyr_RadDwFlux )
      xyr_RadSUwFlux = xyr_RadSUwFlux + xyr_RadUwFlux
      xyr_RadSDwFlux = xyr_RadSDwFlux + xyr_RadDwFlux
    end do
    ! * 1000 to 14286 cm-1 (0.70-10 micron): 
    !
    SolarFluxTOA = SolarConst / DistFromStarScld**2 * DiurnalMeanFactor
    TotUVVISFracSolarFlux = 0.0_DP
    do l = 1, nbands1
      ! UV and Visible optical properties and solar flux
      !
      call RadCL1996UVVISParams( l, UVVISFracSolarFlux, UVVISO3AbsCoef, UVVISRayScatCoef )
      TotUVVISFracSolarFlux = TotUVVISFracSolarFlux + UVVISFracSolarFlux
    end do
    if ( TotUVVISFracSolarFlux > 1.0_DP ) then
      call MessageNotify( 'E', module_name, 'TotUVVISFracSolarFlux, %f, is greater than 1.', d = (/TotUVVISFracSolarFlux/) )
    end if
    SolarFluxTOA = SolarFluxTOA * ( 1.0_DP - TotUVVISFracSolarFlux )
    n = 1
    xyz_DelOptDepH2OVap = a_AbsCoefH2OVap(n) * ( xyz_Press / a_RefPressH2OVap(n) )**a_PressScaleIndH2OVap(n) * xyz_DelAtmMass * xyz_QH2OVap
    !
    xyr_TotOptDep = 0.0_DP
    do k = kmax, 1, -1
      xyr_TotOptDep(:,:,k-1) = xyr_TotOptDep(:,:,k) + xyz_DelOptDepH2OVap(:,:,k)
    end do
!!$    xyz_SSA = 1.0_DP - 1.0e-100_DP
!!$    xyz_AF  = 0.0_DP
!!$
!!$    call RadRTETwoStreamAppSW(        &
!!$      & xyz_SSA, xyz_AF,              & ! (in)
!!$      & xyr_TotOptDep,                & ! (in)
!!$      & xy_SurfAlbedo,                & ! (in)
!!$      & SolarFluxTOA, xy_InAngle,     & ! (in)
!!$      & xyr_RadUwFlux, xyr_RadDwFlux  & ! (out)
!!$      & )
    call RadRTETwoStreamAppHomogAtm( xy_SurfAlbedo, SolarFluxTOA, xy_InAngle, 1.0_DP - 1.0e-10_DP, 0.0_DP, xyr_TotOptDep, xyr_RadUwFlux, xyr_RadDwFlux )
    xyr_RadSUwFlux = xyr_RadSUwFlux + xyr_RadUwFlux
    xyr_RadSDwFlux = xyr_RadSDwFlux + xyr_RadDwFlux
    ! 計算時間計測一時停止
    ! Pause measurement of computation time
    !
    call TimesetClockStop( module_name )
  end subroutine RadSimpleSWV20Flux
          | Subroutine : | 
This procedure input/output NAMELIST#rad_simple_SW_V2_0_nml .
  subroutine RadSimpleSWV20Init
    ! ファイル入出力補助
    ! File I/O support
    !
    use dc_iounit, only: FileOpen
    ! NAMELIST ファイル入力に関するユーティリティ
    ! Utilities for NAMELIST file input
    !
    use namelist_util, only: namelist_filename, NmlutilMsg, NmlutilAryValid
    ! 太陽放射フラックスの設定
    ! Set solar constant
    !
    use set_solarconst, only : SetSolarConstInit
    ! 短波入射 (太陽入射)
    ! Short wave (insolation) incoming
    !
    use rad_short_income, only : RadShortIncomeInit
    !
    ! Solve radiative transfer equation in two stream approximation
    !
    use rad_rte_two_stream_app, only: RadRTETwoStreamAppInit
    ! Chou and Lee (1996) による短波放射モデル
    ! Short wave radiation model described by Chou and Lee (1996)
    !
    use rad_CL1996, only : RadCL1996Init
    ! 宣言文 ; Declaration statements
    !
    real(DP) :: WNBnds             (0:MaxNmlArySize)
    real(DP) :: AbsCoefH2OVap      (1:MaxNmlArySize)
    real(DP) :: PressScaleIndH2OVap(1:MaxNmlArySize)
    real(DP) :: RefPressH2OVap     (1:MaxNmlArySize)
    integer:: unit_nml        ! NAMELIST ファイルオープン用装置番号.
                              ! Unit number for NAMELIST file open
    integer:: iostat_nml      ! NAMELIST 読み込み時の IOSTAT.
                              ! IOSTAT of NAMELIST read
    integer:: n
    ! NAMELIST 変数群
    ! NAMELIST group name
    !
    namelist /rad_simple_SW_V2_0_nml/ FlagRayleighScattering, RayleighScatCoefRatio, AbsCoefH2OVap, PressScaleIndH2OVap, RefPressH2OVap, DiffFact
          !
          ! デフォルト値については初期化手続 "rad_simple_SW_V2_0#RadSimpleSWV20Init"
          ! のソースコードを参照のこと.
          !
          ! Refer to source codes in the initialization procedure
          ! "rad_simple_SW_V2_6#RadSimpleSWEV20Init" for the default values.
          !
    if ( rad_simple_SW_V2_0_inited ) return
    ! デフォルト値の設定
    ! Default values settings
    !
    FlagRayleighScattering = .true.
    RayleighScatCoefRatio  = 1.0_DP
    nbmax = 1
    WNBnds                 = -999.9_DP
    AbsCoefH2OVap          = -999.9_DP
    PressScaleIndH2OVap    = -999.9_DP
    RefPressH2OVap         = -999.9_DP
    AbsCoefH2OVap      (1:nbmax) = (/ 0.0d0  /)
    PressScaleIndH2OVap(1:nbmax) = (/ 0.0d0  /)
    RefPressH2OVap     (1:nbmax) = (/ 1.0d5  /)
    DiffFact        = 1.66_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 = rad_simple_SW_V2_0_nml, iostat = iostat_nml )                  ! (out)
      close( unit_nml )
      call NmlutilMsg( iostat_nml, module_name ) ! (in)
    end if
    a_AbsCoefH2OVap       = AbsCoefH2OVap
    a_PressScaleIndH2OVap = PressScaleIndH2OVap
    a_RefPressH2OVap      = RefPressH2OVap
    ! Initialization of modules used in this module
    !
    ! 太陽放射フラックスの設定
    ! Set solar constant
    !
    call SetSolarConstInit
    ! 短波入射 (太陽入射)
    ! Short wave (insolation) incoming
    !
    call RadShortIncomeInit
    !
    ! Solve radiative transfer equation in two stream approximation
    !
    call RadRTETwoStreamAppInit
    ! Chou and Lee (1996) による短波放射モデル
    ! Short wave radiation model described by Chou and Lee (1996)
    !
    call RadCL1996Init
    ! 印字 ; Print
    !
    call MessageNotify( 'M', module_name, '----- Initialization Messages -----' )
    call MessageNotify( 'M', module_name, 'FlagRayleighScattering = %b', l = (/ FlagRayleighScattering /) )
    call MessageNotify( 'M', module_name, 'RayleighScatCoefRatio  = %f', d = (/ RayleighScatCoefRatio /) )
    do n = 1, nbmax
!!$      call MessageNotify( 'M', module_name, '  %d : %f %f %f %f %f %f %f %f',       &
      call MessageNotify( 'M', module_name, '  %d : %f %f %f %f %f', i = (/ n /), d = (/ a_WNBnds(n-1)*1.0e-2, a_WNBnds(n)*1.0e-2, a_AbsCoefH2OVap(n), a_PressScaleIndH2OVap(n), a_RefPressH2OVap(n) /) )
    end do
    call MessageNotify( 'M', module_name, 'DiffFact        = %f', d = (/ DiffFact /) )
    !
    call MessageNotify( 'M', module_name, '-- version = %c', c1 = trim(version) )
    rad_simple_SW_V2_0_inited = .true.
  end subroutine RadSimpleSWV20Init
          | Variable : | |||
| RayleighScatCoefRatio : | real(DP), save
  | 
| Variable : | |||
| a_AbsCoefH2OVap(1:MaxNmlArySize) : | real(DP), save
  | 
| Variable : | |||
| a_PressScaleIndH2OVap(1:MaxNmlArySize) : | real(DP), save
  | 
| Variable : | |||
| a_RefPressH2OVap(1:MaxNmlArySize) : | real(DP), save
  | 
| Constant : | |||
| module_name = ‘rad_simple_SW_V2_0‘ : | character(*), parameter
  | 
| Variable : | |||
| rad_simple_SW_V2_0_inited = .false. : | logical, save
  | 
| Constant : | |||
| version = ’$Name: $’ // ’$Id: rad_Earth_SW_V2_6.f90,v 1.1 2015/01/29 12:16:19 yot Exp $’ : | character(*), parameter
  |