| Class | rad_SL09 | 
| In: | 
                
                radiation/rad_SL09.f90
                
         | 
Note that Japanese and English are described in parallel.
This is a radiation model described by Schneider and Liu (2009).
Schneider, T. and J. Liu, Formation of jets and equatorial superrotation on Jupiter, J. Atmos. Sci., 69, 579, 2009.
| !$ ! 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) | 
| Subroutine : | |
| xyr_Press( 0:imax-1, 1:jmax, 0:kmax ) : | real(DP), intent(in ) | 
| xyz_Temp( 0:imax-1, 1:jmax, 1:kmax ) : | real(DP), intent(in ) | 
| xyr_RadSFlux( 0:imax-1, 1:jmax, 0:kmax ) : | real(DP), intent(out) | 
| xyr_RadLFlux(0:imax-1, 1:jmax, 0:kmax) : | real(DP), intent(out) | 
| xyra_DelRadLFlux(0:imax-1, 1:jmax, 0:kmax, 0:1) : | real(DP), intent(out) | 
  subroutine RadSL09Flux( xyr_Press, xyz_Temp, xyr_RadSFlux, xyr_RadLFlux, xyra_DelRadLFlux )
    ! USE statements
    !
    !
    ! Physical constants settings
    !
    use constants, only: PI       ! $ \pi $ .
                                  ! Circular constant
    ! 座標データ設定
    ! Axes data settings
    !
    use axesset, only : y_Lat
    !
    ! Solve radiative transfer equation in two stream approximation
    !
    use rad_two_stream_app, only: RadTwoStreamApp, IDSchemeShortWave, RadTwoStreamAppHomogAtm
    real(DP), intent(in ) :: xyr_Press    ( 0:imax-1, 1:jmax, 0:kmax )
    real(DP), intent(in ) :: xyz_Temp     ( 0:imax-1, 1:jmax, 1:kmax )
    real(DP), intent(out) :: xyr_RadSFlux ( 0:imax-1, 1:jmax, 0:kmax )
    real(DP), intent(out) :: xyr_RadLFlux    (0:imax-1, 1:jmax, 0:kmax)
    real(DP), intent(out) :: xyra_DelRadLFlux(0:imax-1, 1:jmax, 0:kmax, 0:1)
    ! Work variables
    !
    real(DP) :: SolarFluxTOA
!!$    real(DP) :: QeRatio
!!$    real(DP) :: xyz_SSA      (0:imax-1, 1:jmax, 1:kmax)
!!$    real(DP) :: xyz_AF       (0:imax-1, 1:jmax, 1:kmax)
    real(DP) :: xy_SurfAlbedo(0:imax-1, 1:jmax)
    real(DP) :: xy_InAngle   (0:imax-1, 1:jmax)
    real(DP) :: xy_CosSZA    (0:imax-1, 1:jmax)
    real(DP) :: xyr_OptDep   (0:imax-1, 1:jmax, 0:kmax)
    real(DP) :: SSA
    real(DP) :: AF
    integer  :: i
    integer  :: j
    integer  :: k
    ! 初期化
    ! Initialization
    !
    if ( .not. rad_SL09_inited ) call RadSL09Init
    ! Short wave radiation
    !
    xyr_OptDep = SWOptDepAtRefPress * ( xyr_Press / SWRefPress )**SWOrd
    SolarFluxTOA  = SolarConst / PI
    SSA = 0.8_DP
    AF  = 0.204_DP
    !   Af = 0 may be much better than 0.204 when Eddington approximation is used.
!!$    AF         = 0.0_DP
    do j = 1, jmax
      xy_CosSZA(:,j) = cos( y_Lat(j) )
    end do
    do j = 1, jmax
      do i = 0, imax-1
        if ( xy_CosSZA(i,j) > 0.0_DP ) then
          xy_InAngle(i,j) = 1.0_DP / xy_CosSZA(i,j)
        else
          xy_InAngle(i,j) = 0.0_DP
        end if
      end do
    end do
    !   Unused variable but this is required as an argument
    !
    xy_SurfAlbedo = 1.0d100
    call RadTwoStreamAppHomogAtm( xy_SurfAlbedo, SolarFluxTOA, xy_InAngle, SSA, AF, xyr_OptDep, xyr_RadSFlux, FlagSemiInfAtm = .true., FlagSL09 = .true. )
    ! Long wave radiation
    !
    !   Although the surface temperature and surface emissivity are set, but are not used.
    !
    xyr_OptDep = LWOptDepAtRefPress * ( xyr_Press / LWRefPress )**LWOrd
!!$    call RadiationRTEQNonScat(                    &
!!$      & xyz_Temp, xy_SurfTemp, xy_SurfEmis, xyr_OptDep, & ! (in )
!!$      & xyr_RadLFlux, xyra_DelRadLFlux,                 & ! (out)
!!$      & xy_SurfUpRadLFluxBase = xyr_RadSFlux(:,:,0)     & ! (in ) optional
!!$      & )
    call RadSL09LWFlux( xyz_Temp, xyr_OptDep, xyr_RadSFlux(:,:,0), xyr_RadLFlux, xyra_DelRadLFlux )
  end subroutine RadSL09Flux
          | Variable : | |||
| rad_SL09_inited = .false. : | logical, save, public
  | 
| Subroutine : | 
This procedure input/output NAMELIST#rad_SL09_nml .
  subroutine RadSL09Init
    ! ファイル入出力補助
    ! File I/O support
    !
    use dc_iounit, only: FileOpen
    ! NAMELIST ファイル入力に関するユーティリティ
    ! Utilities for NAMELIST file input
    !
    use namelist_util, only: namelist_filename, NmlutilMsg, NmlutilAryValid
    ! メッセージ出力
    ! Message output
    !
    use dc_message, only: MessageNotify
    ! 宣言文 ; 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_SL09_nml/ SWOptDepAtRefPress, SWRefPress, SWOrd, LWOptDepAtRefPress, LWRefPress, LWOrd, SolarConst
          !
          ! デフォルト値については初期化手続 "rad_SL09#RadSL09Init"
          ! のソースコードを参照のこと.
          !
          ! Refer to source codes in the initialization procedure
          ! "rad_SL09#RadSL09Init" for the default values.
          !
    ! デフォルト値の設定
    ! Default values settings
    !
    SWOptDepAtRefPress =  3.0_DP
    SWRefPress         =  3.0d5
    SWOrd              =  1.0_DP
    LWOptDepAtRefPress = 80.0_DP
    LWRefPress         =  3.0d5
    LWOrd              =  2.0_DP
    SolarConst         = 50.7_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_SL09_nml, iostat = iostat_nml )             ! (out)
      close( unit_nml )
      call NmlutilMsg( iostat_nml, module_name ) ! (in)
    end if
    ! 印字 ; Print
    !
    call MessageNotify( 'M', module_name, '----- Initialization Messages -----' )
    call MessageNotify( 'M', module_name, 'SWOptDepAtRefPress = %f', d = (/ SWOptDepAtRefPress /) )
    call MessageNotify( 'M', module_name, 'SWRefPress         = %f', d = (/ SWRefPress /) )
    call MessageNotify( 'M', module_name, 'SWOrd              = %f', d = (/ SWOrd /) )
    call MessageNotify( 'M', module_name, 'LWOptDepAtRefPress = %f', d = (/ LWOptDepAtRefPress /) )
    call MessageNotify( 'M', module_name, 'LWRefPress         = %f', d = (/ LWRefPress /) )
    call MessageNotify( 'M', module_name, 'LWOrd              = %f', d = (/ LWOrd /) )
    call MessageNotify( 'M', module_name, 'SolarConst         = %f', d = (/ SolarConst /) )
    call MessageNotify( 'M', module_name, '-- version = %c', c1 = trim(version) )
    rad_SL09_inited = .true.
  end subroutine RadSL09Init
          | Subroutine : | |||
| xyz_Temp(0:imax-1, 1:jmax, 1:kmax) : | real(DP), intent(in )
  | ||
| xyr_OptDep(0:imax-1, 1:jmax, 0:kmax) : | real(DP), intent(in )
  | ||
| xy_SurfUpRadLFluxBase(0:imax-1, 1:jmax) : | real(DP), intent(in ) | ||
| xyr_RadLFlux(0:imax-1, 1:jmax, 0:kmax) : | real(DP), intent(out)
  | ||
| xyra_DelRadLFlux(0:imax-1, 1:jmax, 0:kmax, 0:1) : | real(DP), intent(out)
  | 
散乱なしの場合の放射伝達方程式の計算
Integrate radiative transfer equation without scattering
  subroutine RadSL09LWFlux( xyz_Temp, xyr_OptDep, xy_SurfUpRadLFluxBase, xyr_RadLFlux, xyra_DelRadLFlux )
    !
    ! 散乱なしの場合の放射伝達方程式の計算
    !
    ! Integrate radiative transfer equation without scattering
    !
    ! モジュール引用 ; USE statements
    !
    ! 物理定数設定
    ! Physical constants settings
    !
    use constants, only: PI, StB
                              ! $ \sigma_{SB} $ .
                              ! ステファンボルツマン定数.
                              ! Stefan-Boltzmann constant
    ! 宣言文 ; Declaration statements
    !
    real(DP), intent(in ) :: xyz_Temp    (0:imax-1, 1:jmax, 1:kmax)
                              ! $ T $ .     温度. Temperature
    real(DP), intent(in ) :: xyr_OptDep  (0:imax-1, 1:jmax, 0:kmax)
                              ! Optical depth
    real(DP), intent(in ) :: xy_SurfUpRadLFluxBase(0:imax-1, 1:jmax)
    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), parameter :: DiffFact = 1.66_DP
    real(DP):: xyr_RadLDoFlux (0:imax-1, 1:jmax, 0:kmax)
    real(DP):: xyr_RadLUpFlux (0:imax-1, 1:jmax, 0:kmax)
    real(DP):: xyz_DelTrans (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)
                              ! Integrated Planck function
    real(DP):: xyz_IntDPFDT     (0:imax-1, 1:jmax, 1:kmax)
                              ! Integrated temperature derivative of Planck function
    real(DP):: xy_SurfUpRadLFlux(0:imax-1, 1:jmax)
    integer:: k, kk           ! 鉛直方向に回る DO ループ用作業変数
                              ! Work variables for DO loop in vertical direction
    ! 実行文 ; Executable statement
    !
    ! 初期化
    ! Initialization
    !
    if ( .not. rad_SL09_inited ) call RadSL09Init
    ! Case for grey atmosphere
    !
    xyz_IntPF       = StB * xyz_Temp**4
    xyz_IntDPFDT    = 4.0_DP * StB * xyz_Temp**3
    ! 透過関数計算
    ! Calculate transmission functions
    !
    do k = 1, kmax
      xyz_DelTrans(:,:,k) = exp( - DiffFact * ( xyr_OptDep(:,:,k-1) - xyr_OptDep(:,:,k) ) )
    end do
    !
    do k = 0, kmax
      do kk = k, k
        xyrr_Trans(:,:,k,kk) = 1.0_DP
      end do
      do kk = k+1, kmax
        xyrr_Trans(:,:,k,kk) = xyrr_Trans(:,:,k,kk-1) * xyz_DelTrans(:,:,kk)
      end do
    end do
    do k = 0, kmax
      do kk = 0, k-1
        xyrr_Trans(:,:,k,kk) = xyrr_Trans(:,:,kk,k)
      end do
    end do
    ! 放射フラックス計算
    ! Calculate radiation flux
    !
    !   Initialization
    !
    xyr_RadLDoFlux = 0.0_DP
    xyr_RadLUpFlux = 0.0_DP
    !
    !   Downward flux
    !
    do k = kmax, 0, -1
      do kk = kmax, k+1, -1
        xyr_RadLDoFlux(:,:,k) = xyr_RadLDoFlux(:,:,k) + xyz_IntPF(:,:,kk) * ( xyrr_Trans(:,:,k,kk-1) - xyrr_Trans(:,:,k,kk) )
      end do
    end do
    !
    !   Upward flux
    !
    !     Set upward flux
    !
    xy_SurfUpRadLFlux = xyr_RadLDoFlux(:,:,0) - xy_SurfUpRadLFluxBase
    !
    do k = 0, kmax
      xyr_RadLUpFlux(:,:,k) = xy_SurfUpRadLFlux * xyrr_Trans(:,:,k,0)
      do kk = 1, k
        xyr_RadLUpFlux(:,:,k) = xyr_RadLUpFlux(:,:,k) - xyz_IntPF(:,:,kk) * ( xyrr_Trans(:,:,k,kk-1) - xyrr_Trans(:,:,k,kk) )
      end do
    end do
    xyr_RadLFlux = xyr_RadLUpFlux - xyr_RadLDoFlux
    ! 放射フラックスの変化率の計算
    ! Calculate rate of change of radiative flux
    !
    do k = 0, kmax
      xyra_DelRadLFlux(:,:,k,0) = 0.0_DP
      xyra_DelRadLFlux(:,:,k,1) = - xyz_IntDPFDT(:,:,1) * ( xyrr_Trans(:,:,k,0) - xyrr_Trans(:,:,k,1) )
    end do
  end subroutine RadSL09LWFlux
          | Subroutine : | |
| SolarFluxTOA : | real(DP), intent(in ) | 
| xy_CosSZA(0:imax-1, 1:jmax) : | real(DP), intent(in ) | 
| SSA : | real(DP), intent(in ) | 
| AF : | real(DP), intent(in ) | 
| xyr_OptDep( 0:imax-1, 1:jmax, 0:kmax ) : | real(DP), intent(in ) | 
| xyr_RadSFlux( 0:imax-1, 1:jmax, 0:kmax ) : | real(DP), intent(out) | 
  subroutine RadSL09SWFlux( SolarFluxTOA, xy_CosSZA, SSA, AF, xyr_OptDep, xyr_RadSFlux )
    real(DP), intent(in ) :: SolarFluxTOA
    real(DP), intent(in ) :: xy_CosSZA(0:imax-1, 1:jmax)
    real(DP), intent(in ) :: SSA
    real(DP), intent(in ) :: AF
    real(DP), intent(in ) :: xyr_OptDep   ( 0:imax-1, 1:jmax, 0:kmax )
    real(DP), intent(out) :: xyr_RadSFlux ( 0:imax-1, 1:jmax, 0:kmax )
    ! Work variables
    !
    real(DP) :: BondAlbedo
    real(DP) :: Gamma
    integer  :: j, k
    BondAlbedo = ( sqrt( 1.0_DP - SSA * AF ) - sqrt( 1.0_DP - SSA ) ) / ( sqrt( 1.0_DP - SSA * AF ) + sqrt( 1.0_DP - SSA ) )
    Gamma = 2.0_DP * sqrt( 1.0_DP - SSA ) * sqrt( 1.0_DP - SSA * AF )
    do k = 0, kmax
      do j = 1, jmax
        xyr_RadSFlux(:,j,k) = - SolarFluxTOA * xy_CosSZA(:,j) * ( 1.0_DP - BondAlbedo ) * exp( -Gamma * xyr_OptDep(:,j,k) )
      end do
    end do
  end subroutine RadSL09SWFlux
          | Constant : | |||
| version = ’$Name: dcpam5-20110327 $’ // ’$Id: rad_SL09.f90,v 1.1 2011-03-27 02:20:12 yot Exp $’ : | character(*), parameter
  |