| Class | radiation_dcpam_M_V1 | 
| In: | 
                
                radiation/radiation_dcpam_M_V1.f90
                
         | 
Note that Japanese and English are described in parallel.
地球大気向け放射モデル.
This is a radiation model for the Mars’ atmospehre.
Radiation in the wavenumber range from 40 to 2600 cm-1 is calculated in the routine of long wave radiation. Radiation in the wavenumber range from 2600 to 66667 cm-1 (0.15 to 3.85 micron) is calculated in the routine of shortwave radiation.
| !$ ! 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 : | |
| xy_SurfAlbedo( 0:imax-1, 1:jmax ) : | real(DP), intent(in ) | 
| xyz_Press( 0:imax-1, 1:jmax, 1:kmax ) : | 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(in ) | 
| xy_SurfTemp(0:imax-1, 1:jmax) : | 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 RadiationDcpamMV1Flux( xy_SurfAlbedo, xyz_Press, xyr_Press, xyz_Temp, xy_SurfTemp, xyr_RadSFlux, xyr_RadLFlux, xyra_DelRadLFlux )
    ! USE statements
    !
    ! 物理定数設定
    ! Physical constants settings
    !
    use constants, only : PI, GasRDry
    ! 時刻管理
    ! Time control
    !
    use timeset, only: TimeN, DelTime, TimesetClockStart, TimesetClockStop
    ! 短波入射 (太陽入射)
    ! Short wave (insolation) incoming
    !
    use radiation_short_income, only : ShortIncoming
    use radiation_dcpam_M_15m, only : rad15m_main
    ! 放射関連ルーチン
    ! Routines for radiation calculation
    !
    use radiation_utils, only : RadiationRTEQNonScatWrapper
    use radiation_two_stream_app, only : RadiationTwoStreamAppHomogAtm
    real(DP), intent(in ) :: xy_SurfAlbedo( 0:imax-1, 1:jmax )
    real(DP), intent(in ) :: xyz_Press    ( 0:imax-1, 1:jmax, 1:kmax )
    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(in ) :: xy_SurfTemp     (0:imax-1, 1:jmax)
    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)
    real(DP), parameter :: DiffFact = 1.66_DP
    real(DP) :: xyz_Rho             (0:imax-1, 1:jmax, 1:kmax)
    real(DP) :: xyr_DOD067  (0:imax-1, 1:jmax, 0:kmax)
    real(DP) :: xyr_DOD     (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)
    real(DP) :: QeRat
    real(DP) :: SSA
    real(DP) :: AF
    real(DP) :: xyr_RadLFluxComp    (0:imax-1, 1:jmax, 0:kmax)
    real(DP) :: xyra_DelRadLFluxComp(0:imax-1, 1:jmax, 0:kmax, 0:1)
    real(DP) :: xy_SurfEmis(0:imax-1, 1:jmax)
    real(DP) :: xy_InAngle (0:imax-1, 1:jmax)
    real(DP) :: DistFromStarScld
    real(DP) :: SolarFluxTOA
    real(DP)           :: WNs
    real(DP)           :: WNe
    integer, parameter :: NumGaussNode = 5
    integer :: k
    integer :: kk
    ! 初期化
    ! Initialization
    !
    if ( .not. radiation_dcpam_M_V1_inited ) call RadiationDcpamMV1Init
    xyz_Rho = xyz_Press / ( GasRDry * xyz_Temp )
    !  Dust optical depth
    !
    call SetDOD067( DOD067, xyr_Press, xyz_Press, xyr_DOD067 )
    ! 短波放射
    ! Short wave radiation
    !
    call ShortIncoming( xy_InAngle         = xy_InAngle, DistFromStarScld   = DistFromStarScld )
    SolarFluxTOA = SolarConst / DistFromStarScld**2
!!$    QeRat   = 0.9837_DP    !   Ockert-Bell et al. (1997)
!!$    xyz_SSA = 0.86_DP
!!$    xyz_AF  = 0.68_DP
!!$    QeRat   = 1.0_DP       !   Clancy and Lee (1991)
    SSA = 0.92_DP
    AF  = 0.55_DP
    call RadiationTwoStreamAppHomogAtm( xy_SurfAlbedo, SolarFluxTOA, xy_InAngle, SSA, AF, xyr_DOD067, xyr_RadSFlux )
    ! 長波放射
    ! Long wave radiation
    !
    xyr_RadLFlux     = 0.0_DP
    xyra_DelRadLFlux = 0.0_DP
    !  Surface emissivity
    !
    xy_SurfEmis = 1.0_DP
    !    Flux from 0 to 500 cm-1
    !
    WNs     =   0.0d2
    WNe     = 500.0d2
    QeRat   = 0.17_DP                       ! Wavenumber averaged extinction coefficient
                                            ! ssa = 0.35d0
                                            ! af  = 0.36d0
    xyr_DOD = QeRat * xyr_DOD067
    do k = 1, kmax
      xyz_DelTrans(:,:,k) = exp( - DiffFact * ( xyr_DOD(:,:,k-1) - xyr_DOD(:,:,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
    call RadiationRTEQNonScatWrapper( xyz_Temp, xy_SurfTemp, xy_SurfEmis, xyrr_Trans, xyr_RadLFluxComp, xyra_DelRadLFluxComp, WNs, WNe, NumGaussNode )
    xyr_RadLFlux     = xyr_RadLFlux     + xyr_RadLFluxComp
    xyra_DelRadLFlux = xyra_DelRadLFlux + xyra_DelRadLFluxComp
    !    Flux from 500 to 850 cm-1
    !
    QeRat = 0.25_DP                      ! Wavenumber averaged extinction coefficient
    SSA   = 0.45_DP                      ! Wavenumber averaged single scattering albedo
    call rad15m_main( TimeN, DelTime, xyz_Temp, xyr_Press, xyz_Press, xy_SurfTemp, xyz_Rho, xyr_DOD067, QeRat, SSA, xy_SurfEmis, xyr_RadLFluxComp, xyra_DelRadLFluxComp )
    xyr_RadLFlux     = xyr_RadLFlux     + xyr_RadLFluxComp
    xyra_DelRadLFlux = xyra_DelRadLFlux + xyra_DelRadLFluxComp
    !    Flux from 850 to 2000 cm-1
    !
    WNs     =  850.0d2
    WNe     = 2000.0d2
    QeRat   =    0.41_DP                    ! Wavenumber averaged extinction coefficient
                                            ! ssa = 0.55d0
                                            ! af  = 0.55d0
    xyr_DOD = QeRat * xyr_DOD067
    do k = 1, kmax
      xyz_DelTrans(:,:,k) = exp( - DiffFact * ( xyr_DOD(:,:,k-1) - xyr_DOD(:,:,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
    call RadiationRTEQNonScatWrapper( xyz_Temp, xy_SurfTemp, xy_SurfEmis, xyrr_Trans, xyr_RadLFluxComp, xyra_DelRadLFluxComp, WNs, WNe, NumGaussNode )
    xyr_RadLFlux     = xyr_RadLFlux     + xyr_RadLFluxComp
    xyra_DelRadLFlux = xyra_DelRadLFlux + xyra_DelRadLFluxComp
  end subroutine RadiationDcpamMV1Flux
          | Variable : | |||
| radiation_dcpam_M_V1_inited = .false. : | logical, save, public
  | 
| Subroutine : | 
This procedure input/output NAMELIST#radiation_dcpam_M_V1_nml .
  subroutine RadiationDcpamMV1Init
    ! ファイル入出力補助
    ! 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 /radiation_dcpam_M_V1_nml/ SolarConst, DOD067, DustVerDistCoef
          !
          ! デフォルト値については初期化手続 "radiation_dcpam_SWEV1#RadiationDcpamSWEV1Init"
          ! のソースコードを参照のこと.
          !
          ! Refer to source codes in the initialization procedure
          ! "radiation_dcpam_SWEV1#RadiationDcpamSWEV1Init" for the default values.
          !
    ! デフォルト値の設定
    ! Default values settings
    !
    SolarConst      = 1380.0_DP / 1.52_DP**2
    DOD067          = 0.0_DP
    DustVerDistCoef = 0.01_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_dcpam_M_V1_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, 'SolarConst = %f', d = (/ SolarConst /) )
    call MessageNotify( 'M', module_name, 'DOD067     = %f', d = (/ DOD067     /) )
    call MessageNotify( 'M', module_name, '-- version = %c', c1 = trim(version) )
    radiation_dcpam_M_V1_inited = .true.
  end subroutine RadiationDcpamMV1Init
          | Subroutine : | |||
| DOD067 : | 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 )
  | ||
| xyr_DOD067(0:imax-1, 1:jmax, 0:kmax) : | real(DP), intent(out)
  | 
Set dust optical depth at 0.67 micron
  subroutine SetDOD067( DOD067, xyr_Press, xyz_Press, xyr_DOD067 )
    !
    ! 
    !
    ! Set dust optical depth at 0.67 micron
    !
    ! モジュール引用 ; USE statements
    !
    ! 物理定数設定
    ! Physical constants settings
    !
    use constants, only: Grav
    ! 宣言文 ; Declaration statements
    !
    real(DP), intent(in ):: DOD067
                              ! Dust optical depth at 0.67 micron
    real(DP), intent(in ):: xyr_Press    (0:imax-1, 1:jmax, 0:kmax)
                              ! Pressure
    real(DP), intent(in ):: xyz_Press    (0:imax-1, 1:jmax, 1:kmax)
                              ! Pressure
    real(DP), intent(out):: xyr_DOD067   (0:imax-1, 1:jmax, 0:kmax)
                              ! Optical depth
    ! 作業変数
    ! Work variables
    !
    real(DP)            :: xyz_MixRtDust(0:imax-1, 1:jmax, 1:kmax)
    real(DP)            :: xy_DODFac    (0:imax-1, 1:jmax)
    real(DP), parameter :: DustOptDepRefPress  = 610.0_DP
    real(DP), parameter :: DustVerDistRefPress = 610.0_DP
    real(DP)            :: MixRtDust0
    integer:: k               ! 鉛直方向に回る DO ループ用作業変数
                              ! Work variables for DO loop in vertical direction
    ! 実行文 ; Executable statement
    !
    ! 初期化
    ! Initialization
    !
    if ( .not. radiation_dcpam_M_V1_inited ) call RadiationDcpamMV1Init
    MixRtDust0      =   1.0_DP
    xyz_MixRtDust = MixRtDust0 * exp( DustVerDistCoef * ( 1.0_DP - ( DustVerDistRefPress / xyz_Press ) ) )
    k = kmax
    xyr_DOD067(:,:,k) = 0.0_DP
    do k = kmax-1, 0, -1
      xyr_DOD067(:,:,k) = xyr_DOD067(:,:,k+1) + xyz_MixRtDust(:,:,k+1) * ( xyr_Press(:,:,k) - xyr_Press(:,:,k+1) ) / Grav
    end do
    xy_DODFac = DOD067 * xyr_Press(:,:,0) / DustOptDepRefPress / xyr_DOD067(:,:,0)
    do k = 0, kmax
      xyr_DOD067(:,:,k) = xyr_DOD067(:,:,k) * xy_DODFac
    end do
  end subroutine SetDOD067
          | Variable : | |||
| SolarConst : | real(DP), save
  | 
| Constant : | |||
| module_name = ‘radiation_dcpam_M_V1‘ : | character(*), parameter
  | 
| Constant : | |||
| version = ’$Name: dcpam5-20110221-4 $’ // ’$Id: radiation_dcpam_M_V1.f90,v 1.5 2011-02-18 04:38:54 yot Exp $’ : | character(*), parameter
  |