Class rad_SL09
In: radiation/rad_SL09.f90

Schneider and Liu (2009) ���¾å��¢ã����

Radiation model by Schneider and Liu (2009)

Note that Japanese and English are described in parallel.

This is a radiation model described by Schneider and Liu (2009).

References

 Schneider, T. and J. Liu,
   Formation of jets and equatorial superrotation on Jupiter,
   J. Atmos. Sci., 69, 579, 2009.

Procedures List

!$ ! 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)

NAMELIST

NAMELIST#rad_SL09_nml

Methods

Included Modules

dc_types gridset dc_message constants0 axesset rad_rte_two_stream_app dc_iounit namelist_util

Public Instance methods

Subroutine :
xyr_Press(0:imax-1, 1:jmax, 0:kmax) :real(DP), intent(in )
xyz_Press(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(in )
xyz_Temp(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)
xyr_RadLUwFlux(0:imax-1, 1:jmax, 0:kmax) :real(DP), intent(out)
xyr_RadLDwFlux(0:imax-1, 1:jmax, 0:kmax) :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) :real(DP), intent(out)

[Source]

  subroutine RadSL09Flux( xyr_Press, xyz_Press, xyz_Temp, xyr_RadSUwFlux, xyr_RadSDwFlux, xyr_RadLUwFlux, xyr_RadLDwFlux, xyra_DelRadLUwFlux, xyra_DelRadLDwFlux )


    ! USE statements
    !

    ! �¡ã���»ã�¼ã�¸å�ºå��
    ! Message output
    !
    use dc_message, only: MessageNotify

    ! �����»æ�°å­¦å®��°è¨­å®�
    ! Physical and mathematical constants settings
    !
    use constants0, only: PI                    ! $ \pi $ .
                              ! �����.  Circular constant

    ! 座æ����¼ã�¿è¨­å®�
    ! Axes data settings
    !
    use axesset, only : y_Lat

    !
    ! Solve radiative transfer equation in two stream approximation
    !
!!$    use rad_rte_two_stream_app, only: OLD_RadRTETwoStreamAppHomogAtm
    use rad_rte_two_stream_app, only: RadRTETwoStreamAppHomogAtm

    real(DP), intent(in ) :: xyr_Press         (0:imax-1, 1:jmax, 0:kmax)
    real(DP), intent(in ) :: xyz_Press         (0:imax-1, 1:jmax, 1:kmax)
    real(DP), intent(in ) :: xyz_Temp          (0:imax-1, 1:jmax, 1:kmax)
    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), intent(out) :: xyr_RadLUwFlux    (0:imax-1, 1:jmax, 0:kmax)
    real(DP), intent(out) :: xyr_RadLDwFlux    (0:imax-1, 1:jmax, 0:kmax)
    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) :: 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 ) then
      call MessageNotify( 'E', module_name, 'This module has not been initialized.' )
    end if


    ! Short wave radiation
    !
    xyr_OptDep = SWOptDepAtRefPress * ( xyr_Press / SWRefPress )**SWOrd


    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


    if ( FlagGMIns ) then
      SolarFluxTOA = SolarConst / 4.0_DP
      xy_CosSZA    = 1.0_DP
    else
      SolarFluxTOA = SolarConst / PI
      do j = 1, jmax
        xy_CosSZA(:,j) = cos( y_Lat(j) )
      end do
    end if
!!$    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 OLD_RadRTETwoStreamAppHomogAtm(                                   &
!!$      & xy_SurfAlbedo, SolarFluxTOA, xy_InAngle, SSA, AF, xyr_OptDep,  & ! (in )
!!$      & xyr_RadSFlux,                                                  & ! (out)
!!$      & FlagSemiInfAtm = .true., FlagSL09 = .true.                     & ! (in ) optional
!!$      & )
!!$
!!$    call RadRTETwoStreamAppHomogAtm(                               &
!!$      & xy_SurfAlbedo, SolarFluxTOA, xy_InAngle, SSA, AF, xyr_OptDep,  & ! (in )
!!$      & xyr_RadSUwFlux, xyr_RadSDwFlux,                                & ! (out)
!!$      & FlagSemiInfAtm = .true., FlagSL09 = .true                      & ! (in ) optional
!!$      & )

    call RadSL09SWFlux( SolarFluxTOA, xy_CosSZA, SSA, AF, xyr_OptDep, xyr_RadSUwFlux )
    xyr_RadSDwFlux = 0.0_DP



    ! 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,                     & ! (in )
!!$      & xyr_RadSFlux(:,:,0),                      & ! (in )
!!$      & xyr_RadLFlux, xyra_DelRadLFlux            & ! (out)
!!$      & )
    call RadSL09LWFlux( xyr_Press, xyz_Press, xyz_Temp, xyr_OptDep, xyr_RadSUwFlux(:,:,0)-xyr_RadSDwFlux(:,:,0), xyr_RadLUwFlux, xyr_RadLDwFlux, xyra_DelRadLUwFlux, xyra_DelRadLDwFlux )


  end subroutine RadSL09Flux
Subroutine :

This procedure input/output NAMELIST#rad_SL09_nml .

[Source]

  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

    !
    ! Solve radiative transfer equation in two stream approximation
    !
    use rad_rte_two_stream_app, only: RadRTETwoStreamAppInit


    ! 宣�� ; 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/ FlagGMIns, 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.
          !


    if ( rad_SL09_inited ) return


    ! �����������¤ã��¨­å®�
    ! Default values settings
    !
    FlagGMIns          = .false.

    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


    ! Initialization of modules used in this module
    !

    !
    ! Solve radiative transfer equation in two stream approximation
    !
    call RadRTETwoStreamAppInit


    ! �°å� ; Print
    !
    call MessageNotify( 'M', module_name, '----- Initialization Messages -----' )
    call MessageNotify( 'M', module_name, 'FlagGMIns          = %b', l = (/ FlagGMIns /) )
    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
rad_SL09_inited
Variable :
rad_SL09_inited = .false. :logical, save, public
: ����設������. Initialization flag

Private Instance methods

FlagGMIns
Variable :
FlagGMIns :logical, save
LWOptDepAtRefPress
Variable :
LWOptDepAtRefPress :real(DP), save
LWOrd
Variable :
LWOrd :real(DP), save
LWRefPress
Variable :
LWRefPress :real(DP), save
Subroutine :
xyz_Temp(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(in )
: $ T $ . æ¸�º¦. Temperature
xyr_OptDep(0:imax-1, 1:jmax, 0:kmax) :real(DP), intent(in )
: Optical depth
xy_SurfUpRadLFluxBase(0:imax-1, 1:jmax) :real(DP), intent(in )
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

�£ä¹±�������´å�����¾å�ä¼����¹ç�å¼���ç®�

Integrate radiative transfer equation without scattering

[Source]

  subroutine OLD_RadSL09LWFlux( xyz_Temp, xyr_OptDep, xy_SurfUpRadLFluxBase, xyr_RadLFlux, xyra_DelRadLFlux )
    !
    ! �£ä¹±�������´å�����¾å�ä¼����¹ç�å¼���ç®�
    !
    ! Integrate radiative transfer equation without scattering
    !

    ! �¢ã�¸ã�¥ã�¼ã����� ; USE statements
    !

    ! �¡ã���»ã�¼ã�¸å�ºå��
    ! Message output
    !
    use dc_message, only: MessageNotify

    ! �����»æ�°å­¦å®��°è¨­å®�
    ! Physical and mathematical constants settings
    !
    use constants0, 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 ) then
      call MessageNotify( 'E', module_name, 'This module has not been initialized.' )
    end if


    ! 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 OLD_RadSL09LWFlux
Subroutine :
xyr_Press(0:imax-1, 1:jmax, 0:kmax) :real(DP), intent(in )
xyz_Press(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(in )
xyz_Temp(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(in )
: $ T $ . æ¸�º¦. Temperature
xyr_OptDep(0:imax-1, 1:jmax, 0:kmax) :real(DP), intent(in )
: Optical depth
xy_SurfRadSFlux(0:imax-1, 1:jmax) :real(DP), intent(in )
xyr_RadLUwFlux(0:imax-1, 1:jmax, 0:kmax) :real(DP), intent(out)
: �·æ³¢����������. Longwave flux
xyr_RadLDwFlux(0:imax-1, 1:jmax, 0:kmax) :real(DP), intent(out)
: �·æ³¢����������. Longwave flux
xyra_DelRadLUwFlux(0:imax-1, 1:jmax, 0:kmax, 0:1) :real(DP), intent(out)
: �·æ³¢�°è¡¨æ¸�º¦å¤���. Surface temperature tendency with longwave
xyra_DelRadLDwFlux(0:imax-1, 1:jmax, 0:kmax, 0:1) :real(DP), intent(out)
: �·æ³¢�°è¡¨æ¸�º¦å¤���. Surface temperature tendency with longwave

�£ä¹±�������´å�����¾å�ä¼����¹ç�å¼���ç®�

Integrate radiative transfer equation without scattering

[Source]

  subroutine RadSL09LWFlux( xyr_Press, xyz_Press, xyz_Temp, xyr_OptDep, xy_SurfRadSFlux, xyr_RadLUwFlux, xyr_RadLDwFlux, xyra_DelRadLUwFlux, xyra_DelRadLDwFlux )
    !
    ! �£ä¹±�������´å�����¾å�ä¼����¹ç�å¼���ç®�
    !
    ! Integrate radiative transfer equation without scattering
    !

    ! �¢ã�¸ã�¥ã�¼ã����� ; USE statements
    !

    ! �¡ã���»ã�¼ã�¸å�ºå��
    ! Message output
    !
    use dc_message, only: MessageNotify

    ! �����»æ�°å­¦å®��°è¨­å®�
    ! Physical and mathematical constants settings
    !
    use constants0, only: StB
                              ! $ \sigma_{SB} $ .
                              ! �¹ã�����¡ã�³ã���������³å���.
                              ! Stefan-Boltzmann constant

    ! 宣�� ; Declaration statements
    !

    real(DP), intent(in ) :: xyr_Press   (0:imax-1, 1:jmax, 0:kmax)
    real(DP), intent(in ) :: xyz_Press   (0:imax-1, 1:jmax, 1:kmax)
    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_SurfRadSFlux(0:imax-1, 1:jmax)
    real(DP), intent(out) :: xyr_RadLUwFlux (0:imax-1, 1:jmax, 0:kmax)
                              ! �·æ³¢����������. 
                              ! Longwave flux
    real(DP), intent(out) :: xyr_RadLDwFlux (0:imax-1, 1:jmax, 0:kmax)
                              ! �·æ³¢����������. 
                              ! Longwave flux
    real(DP), intent(out) :: xyra_DelRadLUwFlux (0:imax-1, 1:jmax, 0:kmax, 0:1)
                              ! �·æ³¢�°è¡¨æ¸�º¦å¤���. 
                              ! Surface temperature tendency with longwave
    real(DP), intent(out) :: xyra_DelRadLDwFlux (0:imax-1, 1:jmax, 0:kmax, 0:1)
                              ! �·æ³¢�°è¡¨æ¸�º¦å¤���. 
                              ! Surface temperature tendency with longwave


    ! �業��
    ! Work variables
    !
    real(DP):: xyr_Temp       (0:imax-1, 1:jmax, 0:kmax)
    real(DP):: xyr_IntPF      (0:imax-1, 1:jmax, 0:kmax)
                              ! Integrated Planck function
    real(DP):: xyz_DPFDOptDep (0:imax-1, 1:jmax, 1:kmax)
    real(DP):: xyz_DelOptDep     (0:imax-1, 1:jmax, 1:kmax)
    real(DP):: xyz_TransEachLayer(0:imax-1, 1:jmax, 1:kmax)
                              ! ������. 
                              ! Transmission coefficient

    integer:: k               ! ���´æ�¹å�������� DO ���¼ã�����業å���
                              ! Work variables for DO loop in vertical direction

    ! ���� ; Executable statement
    !

    ! ������
    ! Initialization
    !
    if ( .not. rad_SL09_inited ) then
      call MessageNotify( 'E', module_name, 'This module has not been initialized.' )
    end if


    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)


    xyr_IntPF       = StB * xyr_Temp**4
    !
    do k = 1, kmax
      xyz_DelOptDep(:,:,k) = xyr_OptDep(:,:,k-1) - xyr_OptDep(:,:,k)
    end do
    !
    xyz_TransEachLayer = exp( - xyz_DelOptDep )
    !
    do k = 1, kmax
      xyz_DPFDOptDep(:,:,k) = ( xyr_IntPF(:,:,k-1) - xyr_IntPF(:,:,k) ) / max( xyz_DelOptDep(:,:,k), 1.0d-100 )
    end do


    ! �¾å����������¹è�ç®�
    ! Calculate radiation flux
    !

    !
    !   Downward flux
    !
    k = kmax
    xyr_RadLDwFlux(:,:,k) = 0.0_DP
    do k = kmax-1, 0, -1
      xyr_RadLDwFlux(:,:,k) = ( xyr_RadLDwFlux(:,:,k+1) - xyr_IntPF(:,:,k+1) ) * xyz_TransEachLayer(:,:,k+1) + xyr_IntPF(:,:,k) - xyz_DPFDOptDep(:,:,k+1) * ( 1.0_DP - xyz_TransEachLayer(:,:,k+1) )
    end do
    !
    !   Upward flux
    !
    !     Set upward flux
    k = 0
    xyr_RadLUwFlux(:,:,k) = xyr_RadLDwFlux(:,:,0) - xy_SurfRadSFlux
    !
    do k = 1, kmax
      xyr_RadLUwFlux(:,:,k) = ( xyr_RadLUwFlux(:,:,k-1) - xyr_IntPF(:,:,k-1) ) * xyz_TransEachLayer(:,:,k) + xyr_IntPF(:,:,k) + xyz_DPFDOptDep(:,:,k) * ( 1.0_DP - xyz_TransEachLayer(:,:,k) )
    end do


    ! �¾å����������¹ã����������ç®�
    ! Calculate rate of change of radiative flux
    !
    xyra_DelRadLUwFlux = 0.0_DP
    xyra_DelRadLDwFlux = 0.0_DP


  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)

[Source]

  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
SWOptDepAtRefPress
Variable :
SWOptDepAtRefPress :real(DP), save
SWOrd
Variable :
SWOrd :real(DP), save
SWRefPress
Variable :
SWRefPress :real(DP), save
SolarConst
Variable :
SolarConst :real(DP), save
module_name
Constant :
module_name = ‘rad_SL09 :character(*), parameter
: �¢ã�¸ã�¥ã�¼ã������ç§�. Module name
version
Constant :
version = ’$Name: $’ // ’$Id: rad_SL09.f90,v 1.6 2013/05/25 06:49:44 yot Exp $’ :character(*), parameter
: �¢ã�¸ã�¥ã�¼ã�������¼ã�¸ã�§ã�� Module version