Class radiation_dcpam_M_15m
In: radiation/radiation_dcpam_M_15m.f90

Methods

Included Modules

dc_types gridset constants dc_message ckd_module planck_func namelist_util dc_iounit

Public Instance methods

Subroutine :
time :real(DP) , intent(in )
dt :real(DP) , intent(in )
gt(0:imax-1, 1:jmax, 1:kmax) :real(DP) , intent(in )
gph(0:imax-1, 1:jmax, 0:kmax) :real(DP) , intent(in )
gp(0:imax-1, 1:jmax, 1:kmax) :real(DP) , intent(in )
gts(0:imax-1, 1:jmax) :real(DP) , intent(in )
grho(0:imax-1, 1:jmax, 1:kmax) :real(DP) , intent(in )
dod067(0:imax-1, 1:jmax, 0:kmax) :real(DP) , intent(in )
qerat :real(DP) , intent(in )
ssa :real(DP) , intent(in )
emis(0:imax-1, 1:jmax) :real(DP) , intent(in )
gr15mnetflh(0:imax-1, 1:jmax, 0:kmax) :real(DP) , intent(out)
gdr15mnetfldtsh(0:imax-1, 1:jmax, 0:kmax, 0:1) :real(DP) , intent(out)

[Source]

  subroutine rad15m_main( time, dt, gt, gph, gp, gts, grho, dod067, qerat, ssa, emis, gr15mnetflh, gdr15mnetfldtsh )

    real(DP)    , intent(in ) :: time
    real(DP)    , intent(in ) :: dt
    real(DP)    , intent(in ) :: gt             (0:imax-1, 1:jmax, 1:kmax)
    real(DP)    , intent(in ) :: gph            (0:imax-1, 1:jmax, 0:kmax)
    real(DP)    , intent(in ) :: gp             (0:imax-1, 1:jmax, 1:kmax)
    real(DP)    , intent(in ) :: gts            (0:imax-1, 1:jmax)
    real(DP)    , intent(out) :: gr15mnetflh    (0:imax-1, 1:jmax, 0:kmax)
    real(DP)    , intent(out) :: gdr15mnetfldtsh(0:imax-1, 1:jmax, 0:kmax, 0:1)
    real(DP)    , intent(in ) :: grho           (0:imax-1, 1:jmax, 1:kmax)
    real(DP)    , intent(in ) :: dod067         (0:imax-1, 1:jmax, 0:kmax)
    real(DP)    , intent(in ) :: qerat
    real(DP)    , intent(in ) :: ssa
    real(DP)    , intent(in ) :: emis           (0:imax-1, 1:jmax)


    !
    ! local variables
    !
    real(DP) :: gor (0:imax-1, 1:jmax), goru (0:imax-1, 1:jmax), gord (0:imax-1, 1:jmax), gsr (0:imax-1, 1:jmax), gsru (0:imax-1, 1:jmax), gsrd (0:imax-1, 1:jmax)

    ! 実行文 ; Executable statement
    !

    ! 初期化
    ! Initialization
    !
    if ( .not. radiation_dcpam_M_15m_inited ) call rad15m_init

    call rad15m_lowatm_newscheme2006( time, dt, gt, gph, gp, gts, dod067, qerat, ssa, emis, gr15mnetflh, gdr15mnetfldtsh, gor, goru, gord, gsr, gsru, gsrd )

  end subroutine rad15m_main
radiation_dcpam_M_15m_inited
Variable :
radiation_dcpam_M_15m_inited = .false. :logical, save, public
: 初期設定フラグ. Initialization flag.

Private Instance methods

amu
Constant :
amu = 1.6605655d-27 :real(DP), parameter
Subroutine :
gp(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(in )
glnp(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(out)

[Source]

  subroutine calc_lnp( gp, glnp )

    real(DP), intent(in ) :: gp  (0:imax-1, 1:jmax, 1:kmax)
    real(DP), intent(out) :: glnp(0:imax-1, 1:jmax, 1:kmax)


    !
    ! local variables
    !
    integer :: i, j, k


    do k = 1, kmax
      do j = 1, jmax
        do i = 0, imax-1
          glnp(i,j,k) = log( gp(i,j,k) + 1.0d-20 )
        end do
      end do
    end do


  end subroutine calc_lnp
Subroutine :
gph(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(in )
glnph(0:imax-1, 1:jmax, 0:kmax) :real(DP), intent(out)

[Source]

  subroutine calc_lnph( gph, glnph )

    real(DP), intent(in ) :: gph  (0:imax-1, 1:jmax, 1:kmax)
    real(DP), intent(out) :: glnph(0:imax-1, 1:jmax, 0:kmax)


    !
    ! local variables
    !
    integer :: i, j, k


    do k = 0, kmax
      do j = 1, jmax
        do i = 0, imax-1
          glnph(i,j,k) = log( gph(i,j,k) + 1.0d-20 )
        end do
      end do
    end do


  end subroutine calc_lnph
Subroutine :
dlambda :real(DP), intent(in )
emis(0:imax-1, 1:jmax) :real(DP), intent(in )
trans_i2i_toa(0:imax-1, 1:jmax, 0:kmax) :real(DP), intent(in )
: f_{1/2} T_{k+1/2,1/2}
trans_i2i_boa(0:imax-1, 1:jmax, 0:kmax) :real(DP), intent(in )
: f_{km+1/2} T_{k+1/2,km+1/2}
trans_i2i_s(0:imax-1, 1:jmax, 0:kmax) :real(DP), intent(in )
: f_{s} T_{k+1/2,km+1/2}
trans_i2m_lli(0:imax-1, 1:jmax, 0:kmax, 1:kmax) :real(DP), intent(in )
: upper layer interface
trans_i2m_uli(0:imax-1, 1:jmax, 0:kmax, 1:kmax) :real(DP), intent(in )
: lower layer interface
pfh(0:imax-1, 1:jmax, 0:kmax) :real(DP), intent(in )
pfs(0:imax-1, 1:jmax) :real(DP), intent(in )
netflh(0:imax-1, 1:jmax, 0:kmax) :real(DP), intent(out)

[Source]

  subroutine calc_rteq_use_meantrans_arr3d( dlambda, emis, trans_i2i_toa, trans_i2i_boa, trans_i2i_s, trans_i2m_lli, trans_i2m_uli, pfh, pfs, netflh )

    use constants , only : PI

    real(DP), intent(in ) :: dlambda
    real(DP), intent(in ) :: emis(0:imax-1, 1:jmax)
    real(DP), intent(in ) :: trans_i2i_toa(0:imax-1, 1:jmax, 0:kmax)    ! f_{1/2}    T_{k+1/2,1/2}
    real(DP), intent(in ) :: trans_i2i_boa(0:imax-1, 1:jmax, 0:kmax)    ! f_{km+1/2} T_{k+1/2,km+1/2}
    real(DP), intent(in ) :: trans_i2i_s  (0:imax-1, 1:jmax, 0:kmax)    ! f_{s}      T_{k+1/2,km+1/2}
    real(DP), intent(in ) :: trans_i2m_lli(0:imax-1, 1:jmax, 0:kmax, 1:kmax) ! upper layer interface
    real(DP), intent(in ) :: trans_i2m_uli(0:imax-1, 1:jmax, 0:kmax, 1:kmax) ! lower layer interface
    real(DP), intent(in ) :: pfh   (0:imax-1, 1:jmax, 0:kmax)
    real(DP), intent(in ) :: pfs   (0:imax-1, 1:jmax)
    real(DP), intent(out) :: netflh(0:imax-1, 1:jmax, 0:kmax)



    !
    ! local variables
    !
    integer :: i, j, k, k2


    do k = 0, kmax
      do j = 1, jmax
        do i = 0, imax-1
          netflh(i,j,k) = 0.0d0
        end do
      end do
    end do

    do k = 0, kmax

      do j = 1, jmax
        do i = 0, imax-1
          netflh(i,j,k) = netflh(i,j,k) + PI * emis(i,j) * pfs(i,j) * dlambda * trans_i2i_s  (i,j,k) - PI * pfh(i,j,0   ) * dlambda * trans_i2i_boa(i,j,k) + PI * pfh(i,j,kmax) * dlambda * trans_i2i_toa(i,j,k)
        end do
      end do

      do k2 = 1, kmax
        do j = 1, jmax
          do i = 0, imax-1
            netflh(i,j,k) = netflh(i,j,k) - PI * pfh(i,j,k2  ) * dlambda * trans_i2m_uli(i,j,k,k2) + PI * pfh(i,j,k2-1) * dlambda * trans_i2m_lli(i,j,k,k2)
          end do
        end do
      end do

    end do


  end subroutine calc_rteq_use_meantrans_arr3d
Subroutine :
nras :integer , intent(in )
nrps :integer , intent(in )
gph(0:imax-1, 1:jmax, 0:kmax) :real(DP), intent(in )
vmr(0:imax-1, 1:jmax, 1:kmax, 1:nras+nrps) :real(DP), intent(in )
mmmass(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(in )
ac(0:imax-1, 1:jmax, 1:kmax, 1:nras) :real(DP), intent(in )
gdod(0:imax-1, 1:jmax, 0:kmax) :real(DP), intent(in )
trans(0:imax-1, 1:jmax, 0:kmax, 0:kmax) :real(DP), intent(out)

[Source]

  subroutine calc_trans_mp_arr3d( nras, nrps, gph, vmr, mmmass, ac, gdod, trans )

    use constants , only : Grav

    integer , intent(in ) :: nras
    integer , intent(in ) :: nrps
    real(DP), intent(in ) :: gph   (0:imax-1, 1:jmax, 0:kmax)
    real(DP), intent(in ) :: vmr   (0:imax-1, 1:jmax, 1:kmax, 1:nras+nrps)
    real(DP), intent(in ) :: mmmass(0:imax-1, 1:jmax, 1:kmax)
    real(DP), intent(in ) :: ac    (0:imax-1, 1:jmax, 1:kmax, 1:nras)
    real(DP), intent(in ) :: gdod  (0:imax-1, 1:jmax, 0:kmax)
    real(DP), intent(out) :: trans (0:imax-1, 1:jmax, 0:kmax, 0:kmax)


    !
    ! local variables
    !
    real(DP)     :: dopdep(0:imax-1, 1:jmax, 1:kmax)
    real(DP)     :: dtrans(0:imax-1, 1:jmax, 1:kmax)
    real(DP)     :: trans1(0:imax-1, 1:jmax )
    real(DP), parameter :: diffac = 1.66_DP


    integer :: i, j, k, k2, n
    integer :: ks, ke



    do k2 = 0, kmax
      do k = 0, kmax
        do j = 1, jmax
          do i = 0, imax-1
            trans(i,j,k,k2) = 1.0d100
          end do
        end do
      end do
    end do


    do k = 1, kmax
      do j = 1, jmax
        do i = 0, imax-1
          dopdep(i,j,k) = 0.0_DP
        end do
      end do
    end do
    do n = 1, nras
      do k = 1, kmax
        do j = 1, jmax
          do i = 0, imax-1
            dopdep(i,j,k) = dopdep(i,j,k) + ac(i,j,k,n) * vmr(i,j,k,n) / mmmass(i,j,k) * ( gph(i,j,k-1) - gph(i,j,k) ) / Grav
          end do
        end do
      end do
    end do

    !
    ! add dust optical depth
    !
    do k = 1, kmax
      do j = 1, jmax
        do i = 0, imax-1
          dopdep(i,j,k) = dopdep(i,j,k) + gdod(i,j,k-1) - gdod(i,j,k)
        end do
      end do
    end do


    do k = 1, kmax
      do j = 1, jmax
        do i = 0, imax-1
          dtrans(i,j,k) = exp( - dopdep(i,j,k) * diffac )
        end do
      end do
    end do


    !
    ! transmission for "zero thickness" layer ( = 1.0 )
    !
    do ks = 0, kmax
      ke = ks
      do j = 1, jmax
        do i = 0, imax-1
          trans(i,j,ks,ke) = 1.0_DP
        end do
      end do
    end do

    do ks = 0, kmax
      do j = 1, jmax
        do i = 0, imax-1
          trans1(i,j) = 1.0_DP
        end do
      end do
      do ke = ks+1, kmax
        do j = 1, jmax
          do i = 0, imax-1
            trans1(i,j) = trans1(i,j) * dtrans(i,j,ke)
          end do
        end do
        do j = 1, jmax
          do i = 0, imax-1
            trans(i,j,ks,ke) = trans1(i,j)
          end do
        end do
      end do
    end do

    do ks = 0, kmax
      do ke = 0, ks-1
        do j = 1, jmax
          do i = 0, imax-1
            trans(i,j,ks,ke) = trans(i,j,ke,ks)
          end do
        end do
      end do
    end do


  end subroutine calc_trans_mp_arr3d
diffac
Constant :
diffac = 1.0d0 / 1.66d0 :real(DP), parameter
Subroutine :
ks :integer , intent(in )
ke :integer , intent(in )
gt(0:imax-1, 1:jmax, ks:ke) :real(DP), intent(in )
glnp(0:imax-1, 1:jmax, ks:ke) :real(DP), intent(in )
iband :integer , intent(in )
jj(0:imax-1, 1:jmax, ks:ke) :integer , intent(out)
kk(0:imax-1, 1:jmax, ks:ke) :integer , intent(out)

[Source]

  subroutine findindices( ks, ke, gt, glnp, iband, jj, kk )

    use ckd_module, only : ckdp

    integer , intent(in ) :: ks
    integer , intent(in ) :: ke
    real(DP), intent(in ) :: gt(0:imax-1, 1:jmax, ks:ke), glnp(0:imax-1, 1:jmax, ks:ke)
    integer , intent(in ) :: iband
    integer , intent(out) :: jj(0:imax-1, 1:jmax, ks:ke), kk  (0:imax-1, 1:jmax, ks:ke)


    !
    ! local variables
    !
    integer :: i, j, k, l


    do k = ks, ke
      do j = 1, jmax
        do i = 0, imax-1
          kk(i,j,k) = 1
        end do
      end do
    end do

    do l = 1+1, ckdp( iband ) % nt - 1
      do k = ks, ke
        do j = 1, jmax
          do i = 0, imax-1
            if( ckdp( iband ) % t( l ) .le. gt(i,j,k) ) kk(i,j,k) = l
          end do
        end do
      end do
    end do

    do k = ks, ke
      do j = 1, jmax
        do i = 0, imax-1
          jj(i,j,k) = 1
        end do
      end do
    end do
    do l = 1+1, ckdp( iband ) % nlnp - 1
      do k = ks, ke
        do j = 1, jmax
          do i = 0, imax-1
            if( ckdp( iband ) % lnp( l ) .le. glnp(i,j,k) ) jj(i,j,k) = l
          end do
        end do
      end do
    end do


  end subroutine findindices
Subroutine :
gt(0:imax-1, 1:jmax) :real(DP), intent(in )
glnp(0:imax-1, 1:jmax) :real(DP), intent(in )
iband :integer , intent(in )
jj(0:imax-1, 1:jmax) :integer , intent(out)
kk(0:imax-1, 1:jmax) :integer , intent(out)

[Source]

  subroutine findindices2D( gt, glnp, iband, jj, kk )

    real(DP), intent(in ) :: gt(0:imax-1, 1:jmax), glnp(0:imax-1, 1:jmax)
    integer , intent(in ) :: iband
    integer , intent(out) :: jj(0:imax-1, 1:jmax), kk  (0:imax-1, 1:jmax)


    !
    ! local variables
    !
    real(DP) :: gt3d(0:imax-1, 1:jmax, 1:1), glnp3d(0:imax-1, 1:jmax, 1:1)
    integer  :: jj3d(0:imax-1, 1:jmax, 1:1), kk3d  (0:imax-1, 1:jmax, 1:1)


    gt3d  (:,:,1) = gt
    glnp3d(:,:,1) = glnp

    call findindices( 1, 1, gt3d, glnp3d, iband, jj3d, kk3d )

    jj = jj3d(:,:,1)
    kk = kk3d(:,:,1)


  end subroutine findindices2D
Subroutine :
gt(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(in )
glnp(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(in )
iband :integer , intent(in )
jj(0:imax-1, 1:jmax, 1:kmax) :integer , intent(out)
kk(0:imax-1, 1:jmax, 1:kmax) :integer , intent(out)

[Source]

  subroutine findindices3D( gt, glnp, iband, jj, kk )

    real(DP), intent(in ) :: gt(0:imax-1, 1:jmax, 1:kmax), glnp(0:imax-1, 1:jmax, 1:kmax)
    integer , intent(in ) :: iband
    integer , intent(out) :: jj(0:imax-1, 1:jmax, 1:kmax), kk  (0:imax-1, 1:jmax, 1:kmax)


    !
    ! local variables
    !


    call findindices( 1, kmax, gt, glnp, iband, jj, kk )


  end subroutine findindices3D
Subroutine :
gt(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(in )
glnp(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(in )
jj(0:imax-1, 1:jmax, 1:kmax) :integer , intent(in )
kk(0:imax-1, 1:jmax, 1:kmax) :integer , intent(in )
ig :integer , intent(in )
iband :integer , intent(in )
ac(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(out)

[Source]

  subroutine getlnac_givenindices( gt, glnp, jj, kk, ig, iband, ac )

    use ckd_module, only : ckdp

    real(DP), intent(in ) :: gt  (0:imax-1, 1:jmax, 1:kmax)
    real(DP), intent(in ) :: glnp(0:imax-1, 1:jmax, 1:kmax)
    integer , intent(in ) :: jj  (0:imax-1, 1:jmax, 1:kmax)
    integer , intent(in ) :: kk  (0:imax-1, 1:jmax, 1:kmax)
    integer , intent(in ) :: ig
    integer , intent(in ) :: iband
    real(DP), intent(out) :: ac  (0:imax-1, 1:jmax, 1:kmax)


    !
    ! local variables
    !
    real(DP) :: lnac1, lnac2
    integer  :: i, j, k


    do k = 1, kmax
      do j = 1, jmax
        do i = 0, imax-1

          lnac1 = ( ckdp(iband)%lnac( ig, jj(i,j,k)  , kk(i,j,k)+1 ) - ckdp(iband)%lnac( ig, jj(i,j,k)  , kk(i,j,k)   ) ) / ( ckdp(iband)%t( kk(i,j,k)+1 ) - ckdp(iband)%t( kk(i,j,k)   ) ) * ( gt(i,j,k) - ckdp( iband ) % t( kk(i,j,k) ) ) + ckdp(iband)%lnac( ig, jj(i,j,k)  , kk(i,j,k)   )
          lnac2 = ( ckdp(iband)%lnac( ig, jj(i,j,k)+1, kk(i,j,k)+1 ) - ckdp(iband)%lnac( ig, jj(i,j,k)+1, kk(i,j,k)   ) ) / ( ckdp(iband)%t( kk(i,j,k)+1 ) - ckdp(iband)%t( kk(i,j,k)   ) ) * ( gt(i,j,k) - ckdp( iband ) % t( kk(i,j,k) ) ) + ckdp(iband)%lnac( ig, jj(i,j,k)+1, kk(i,j,k)   )

          ac(i,j,k) = ( lnac2 - lnac1 ) / ( ckdp( iband ) % lnp( jj(i,j,k)+1 ) - ckdp( iband ) % lnp( jj(i,j,k)   ) ) * ( glnp(i,j,k) - ckdp( iband ) % lnp( jj(i,j,k) ) ) + lnac1
        end do
      end do
    end do


  end subroutine getlnac_givenindices
Subroutine :
gt(0:imax-1, 1:jmax, 0:kmax) :real(DP), intent(in )
gts(0:imax-1, 1:jmax) :real(DP), intent(in )
iband :integer , intent(in )
pfarr(0:imax-1, 1:jmax, 0:kmax) :real(DP), intent(out)
pfsarr(0:imax-1, 1:jmax) :real(DP), intent(out)

[Source]

  subroutine getpf_arr3d_norat( gt, gts, iband, pfarr, pfsarr )

    use planck_func, only : Integ_PF_GQ_Array3D, Integ_PF_GQ_Array2D

    use ckd_module, only : ckdp

    real(DP), intent(in ) :: gt    (0:imax-1, 1:jmax, 0:kmax)
    real(DP), intent(in ) :: gts   (0:imax-1, 1:jmax)
    integer , intent(in ) :: iband
    real(DP), intent(out) :: pfarr (0:imax-1, 1:jmax, 0:kmax)
    real(DP), intent(out) :: pfsarr(0:imax-1, 1:jmax)


    !
    ! local variables
    !
    integer :: ncp_pfint
    integer :: i, j, k


    ncp_pfint = 5

    call Integ_PF_GQ_Array3D( ckdp(iband)%wnbnds(1), ckdp(iband)%wnbnds(2), ncp_pfint, 0, imax-1, 1, jmax, 0, kmax, gt, pfarr )
    call Integ_PF_GQ_Array2D( ckdp(iband)%wnbnds(1), ckdp(iband)%wnbnds(2), ncp_pfint, 0, imax-1, 1, jmax, gts, pfsarr )

    do k = 0, kmax
      do j = 1, jmax
        do i = 0, imax-1
          pfarr(i,j,k) = pfarr(i,j,k) / ( ckdp(iband)%wnbnds(2) - ckdp(iband)%wnbnds(1) )
        end do
      end do
    end do
    do j = 1, jmax
      do i = 0, imax-1
        pfsarr(i,j) = pfsarr(i,j) / ( ckdp(iband)%wnbnds(2) - ckdp(iband)%wnbnds(1) )
      end do
    end do


  end subroutine getpf_arr3d_norat
Subroutine :
ks :integer , intent(in )
ke :integer , intent(in )
gt(0:imax-1, 1:jmax, ks:ke) :real(DP), intent(in )
glnp(0:imax-1, 1:jmax, ks:ke) :real(DP), intent(in )
jj(0:imax-1, 1:jmax, ks:ke) :integer , intent(in )
kk(0:imax-1, 1:jmax, ks:ke) :integer , intent(in )
ig :integer , intent(in )
iband :integer , intent(in )
pfr(0:imax-1, 1:jmax, ks:ke) :real(DP), intent(out)

[Source]

  subroutine getpfr_givenindices( ks, ke, gt, glnp, jj, kk, ig, iband, pfr )

    use ckd_module, only: ckdp

    integer , intent(in ) :: ks
    integer , intent(in ) :: ke
    real(DP), intent(in ) :: gt  (0:imax-1, 1:jmax, ks:ke)
    real(DP), intent(in ) :: glnp(0:imax-1, 1:jmax, ks:ke)
    integer , intent(in ) :: jj  (0:imax-1, 1:jmax, ks:ke)
    integer , intent(in ) :: kk  (0:imax-1, 1:jmax, ks:ke)
    integer , intent(in ) :: ig, iband
    real(DP), intent(out) :: pfr (0:imax-1, 1:jmax, ks:ke)


    !
    ! local variables
    !
    real(DP) :: pfr1, pfr2
    integer  :: i, j, k, l


    do k = ks, ke
      do j = 1, jmax
        do i = 0, imax-1

          pfr1 = ( ckdp(iband)%pfr( ig, jj(i,j,k)  , kk(i,j,k)+1 ) - ckdp(iband)%pfr( ig, jj(i,j,k)  , kk(i,j,k)   ) ) / ( ckdp(iband)%t( kk(i,j,k)+1 ) - ckdp(iband)%t( kk(i,j,k)   ) ) * ( gt(i,j,k) - ckdp( iband ) % t( kk(i,j,k) ) ) + ckdp(iband)%pfr( ig, jj(i,j,k)  , kk(i,j,k)   )
          pfr2 = ( ckdp(iband)%pfr( ig, jj(i,j,k)+1, kk(i,j,k)+1 ) - ckdp(iband)%pfr( ig, jj(i,j,k)+1, kk(i,j,k)   ) ) / ( ckdp(iband)%t( kk(i,j,k)+1 ) - ckdp(iband)%t( kk(i,j,k)   ) ) * ( gt(i,j,k) - ckdp( iband ) % t( kk(i,j,k) ) ) + ckdp(iband)%pfr( ig, jj(i,j,k)+1, kk(i,j,k)   )


          pfr(i,j,k) = ( pfr2 - pfr1 ) / ( ckdp( iband ) % lnp( jj(i,j,k)+1 ) - ckdp( iband ) % lnp( jj(i,j,k)   ) ) * ( glnp(i,j,k) - ckdp( iband ) % lnp( jj(i,j,k) ) ) + pfr1
        end do
      end do
    end do


  end subroutine getpfr_givenindices
Subroutine :
gt(0:imax-1, 1:jmax) :real(DP), intent(in )
glnp(0:imax-1, 1:jmax) :real(DP), intent(in )
jj(0:imax-1, 1:jmax) :integer , intent(in )
kk(0:imax-1, 1:jmax) :integer , intent(in )
ig :integer , intent(in )
iband :integer , intent(in )
pfr(0:imax-1, 1:jmax) :real(DP), intent(out)

[Source]

  subroutine getpfr_givenindices2D( gt, glnp, jj, kk, ig, iband, pfr )

    real(DP), intent(in ) :: gt  (0:imax-1, 1:jmax)
    real(DP), intent(in ) :: glnp(0:imax-1, 1:jmax)
    integer , intent(in ) :: jj  (0:imax-1, 1:jmax)
    integer , intent(in ) :: kk  (0:imax-1, 1:jmax)
    integer , intent(in ) :: ig, iband
    real(DP), intent(out) :: pfr (0:imax-1, 1:jmax)


    !
    ! local variables
    !
    real(DP) :: gt3d  (0:imax-1, 1:jmax, 1:1)
    real(DP) :: glnp3d(0:imax-1, 1:jmax, 1:1)
    integer  :: jj3d  (0:imax-1, 1:jmax, 1:1)
    integer  :: kk3d  (0:imax-1, 1:jmax, 1:1)
    real(DP) :: pfr3d (0:imax-1, 1:jmax, 1:1)


    gt3d  (:,:,1) = gt
    glnp3d(:,:,1) = glnp
    jj3d  (:,:,1) = jj
    kk3d  (:,:,1) = kk

    call getpfr_givenindices( 1, 1, gt3d, glnp3d, jj3d, kk3d, ig, iband, pfr3d )

    pfr (:,:) = pfr3d (:,:,1)


  end subroutine getpfr_givenindices2D
Subroutine :
gt(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(in )
glnp(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(in )
jj(0:imax-1, 1:jmax, 1:kmax) :integer , intent(in )
kk(0:imax-1, 1:jmax, 1:kmax) :integer , intent(in )
ig :integer , intent(in )
iband :integer , intent(in )
pfr(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(out)

[Source]

  subroutine getpfr_givenindices3D( gt, glnp, jj, kk, ig, iband, pfr )

    use ckd_module, only: ckdp

    real(DP), intent(in ) :: gt  (0:imax-1, 1:jmax, 1:kmax)
    real(DP), intent(in ) :: glnp(0:imax-1, 1:jmax, 1:kmax)
    integer , intent(in ) :: jj  (0:imax-1, 1:jmax, 1:kmax)
    integer , intent(in ) :: kk  (0:imax-1, 1:jmax, 1:kmax)
    integer , intent(in ) :: ig, iband
    real(DP), intent(out) :: pfr (0:imax-1, 1:jmax, 1:kmax)


    !
    ! local variables
    !

    call getpfr_givenindices( 1, kmax, gt, glnp, jj, kk, ig, iband, pfr )


  end subroutine getpfr_givenindices3D
kg1
Variable :
kg1( kg1n ) :real(DP), save
kg1n
Constant :
kg1n = 16 :integer , parameter
kg2
Variable :
kg2( kg2n ) :real(DP), save
kg2n
Constant :
kg2n = 55 :integer , parameter
kg3
Variable :
kg3( kg3n ) :real(DP), save
kg3n
Constant :
kg3n = 3 :integer , parameter
lc
Constant :
lc = 16 :integer , parameter
lnkg
Variable :
lnkg( kg1n, kg2n, kg3n ) :real(DP), save
Subroutine :
m :integer, intent(in )
ig :integer, intent(out)
iband :integer, intent(out)

[Source]

  subroutine m2ckdpindices( m, ig, iband )

    use ckd_module, only : ckdp, nband

    integer, intent(in ) :: m
    integer, intent(out) :: ig
    integer, intent(out) :: iband


    !
    ! local variables
    !
    integer :: num


    ! The comments below will be removed.


    num = 0
    do iband = 1, nband
      if( num + ckdp( iband ) % ng .ge. m ) exit
      num = num + ckdp( iband ) % ng
    end do
    if( iband .gt. nband ) then
      write( 6, * ) 'Unexpected m'
      write( 6, * ) m
      stop
    end if
    ig = m - num
    if( ig .gt. ckdp( iband ) % ng ) then
      write( 6, * ) 'Unexpected ig'
      write( 6, * ) iband, ig
      stop
    end if


  end subroutine m2ckdpindices
module_name
Constant :
module_name = ‘radiation_dcpam_M_15m :character(*), parameter
: モジュールの名称. Module name
nl15fa
Variable :
nl15fa( nl15fn ) :real(DP), save
nl15fn
Constant :
nl15fn = 70 :integer , parameter
nl15sn
Variable :
nl15sn( nl15fn ) :real(DP), save
nlte_refp
Constant :
nlte_refp = 1.0d-2 :real(DP), parameter
nras
Variable :
nras :integer , save
nrps
Variable :
nrps :integer , save
nwnl
Variable :
nwnl :integer , save
nwnsl
Variable :
nwnsl :integer , save
Subroutine :

This procedure input/output NAMELIST#radiation_dcpam_M_15m_nml .

[Source]

  subroutine rad15m_init

    ! モジュール引用 ; 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

    use ckd_module, only : ckd_input, ckdp, nband


    !
    ! local variables
    !
    integer:: unit_nml        ! NAMELIST ファイルオープン用装置番号.
                              ! Unit number for NAMELIST file open
    integer:: iostat_nml      ! NAMELIST 読み込み時の IOSTAT.
                              ! IOSTAT of NAMELIST read

    character(STRING) :: rad15mkg_fn
    character(STRING) :: rad15mnf_fn

    integer           :: m

    namelist /radiation_dcpam_M_15m_nml/ rad15mkg_fn, rad15mint


    ! 実行文 ; Executable statement
    !

    if ( radiation_dcpam_M_15m_inited ) return

    ! デフォルト値の設定
    ! Default values settings
    !

    rad15mkg_fn = "./kg15m"
!!$    rad15mnf_fn = "./nlte15mfactor"
    rad15mint   = 925.0_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_15m_nml, iostat = iostat_nml )         ! (out)
      close( unit_nml )

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


!!$    allocate( rad_gp   ( im, jm,   km ) )
!!$    allocate( rad_gph  ( im, jm, 0:km ) )
!!$    allocate( rad_gt   ( im, jm,   km ) )
!!$    allocate( rad_gts  ( im, jm,   1  ) )
!!$    allocate( rad_gdod ( im, jm, 0:km ) )


    nras = 1
    nrps = 0

!!$    allocate( sgmh_f          ( km*nvr+1 ), &
!!$      &       sgm_f           ( km*nvr   ) )
!!$    allocate( gph_f    ( im, jm, km*nvr+1 ), &
!!$      &       gp_f     ( im, jm, km*nvr   ), &
!!$      &       gth_f    ( im, jm, km*nvr+1 ) )
!!$
!!$    allocate( gvmr_f   ( im, jm, km*nvr  , nras + nrps ) )
!!$    allocate( mmmass_f ( im, jm, km*nvr   ) )
!!$    allocate( ac_f     ( im, jm, km*nvr  , nras        ) )
!!$
!!$    allocate( gdod_f   ( im, jm, km*nvr+1 ) )
!!$
!!$
    allocate( trans  (0:imax-1, 1:jmax, 0:kmax, 0:kmax) )
!!$    allocate( pfh_f    ( im, jm, km*nvr+1 ) )
!!$
!!$    allocate( uwflh_f  ( im, jm, km*nvr+1 ), &
!!$      &       dwflh_f  ( im, jm, km*nvr+1 ) )


    allocate( trans_i2i_toa(0:imax-1, 1:jmax, 0:kmax), trans_i2i_boa(0:imax-1, 1:jmax, 0:kmax), trans_i2i_s  (0:imax-1, 1:jmax, 0:kmax), trans_i2m_lli(0:imax-1, 1:jmax, 0:kmax, 1:kmax), trans_i2m_uli(0:imax-1, 1:jmax, 0:kmax, 1:kmax) )

    trans_i2i_toa(:,:,:)   = 1.0d100
    trans_i2i_boa(:,:,:)   = 1.0d100
    trans_i2i_s  (:,:,:)   = 1.0d100
    trans_i2m_lli(:,:,:,:) = 1.0d100
    trans_i2m_uli(:,:,:,:) = 1.0d100


    !
    ! check
    !
    if( nras .ne. 1 ) then
      write( 6, * ) 'nras is not 1.'
      write( 6, * ) nras
      stop
    end if

    call ckd_input( rad15mkg_fn )

    ! check
    if( nband /= 1 ) then
      write( 6, * ) ' nband is not 1.'
      write( 6, * ) nband
      stop
    end if

    nwnl = 0
    do m = 1, nband
      nwnl = nwnl + ckdp( m ) % ng
    end do



!!$    call increase_vreso_boundary( km, nvr, sgmh, sgmh_f, "log" )
!!$    do k = 1, km * nvr
!!$      sgm_f( k ) = sqrt( sgmh_f( k ) * sgmh_f( k+1 ) )
!!$    end do


!!$    call rad15m_readnlte15mfac( rad15mnf_fn )


    !
    ! This routine must be called after rad15m_readkgtbl.
    !
!!$      call rad15m_rv_read( time )
!!$    call rad15m_rv_read_newscheme2006( time )


    radiation_dcpam_M_15m_inited = .true.


  end subroutine rad15m_init
Subroutine :
time :real(DP) , intent(in )
dt :real(DP) , intent(in )
gt(0:imax-1, 1:jmax, 1:kmax) :real(DP) , intent(in )
gph(0:imax-1, 1:jmax, 0:kmax) :real(DP) , intent(in )
gp(0:imax-1, 1:jmax, 1:kmax) :real(DP) , intent(in )
gts(0:imax-1, 1:jmax) :real(DP) , intent(in )
dod067(0:imax-1, 1:jmax, 0:kmax) :real(DP) , intent(in )
qerat :real(DP) , intent(in )
ssa :real(DP) , intent(in )
emis(0:imax-1, 1:jmax) :real(DP) , intent(in )
gr15mnetflh(0:imax-1, 1:jmax, 0:kmax) :real(DP) , intent(out)
gdr15mnetfldtsh(0:imax-1, 1:jmax, 0:kmax, 0:1) :real(DP) , intent(out)
gor(0:imax-1, 1:jmax) :real(DP) , intent(out)
goru(0:imax-1, 1:jmax) :real(DP) , intent(out)
gord(0:imax-1, 1:jmax) :real(DP) , intent(out)
gsr(0:imax-1, 1:jmax) :real(DP) , intent(out)
gsru(0:imax-1, 1:jmax) :real(DP) , intent(out)
gsrd(0:imax-1, 1:jmax) :real(DP) , intent(out)

[Source]

  subroutine rad15m_lowatm_newscheme2006( time, dt, gt, gph, gp, gts, dod067, qerat, ssa, emis, gr15mnetflh, gdr15mnetfldtsh, gor, goru, gord, gsr, gsru, gsrd )


    use constants , only : Grav, CpDry

    ! メッセージ出力
    ! Message output
    !
    use dc_message, only: MessageNotify

    use ckd_module, only : ckdp


    real(DP)    , intent(in ) :: time
    real(DP)    , intent(in ) :: dt
    real(DP)    , intent(in ) :: gph   (0:imax-1, 1:jmax, 0:kmax)
    real(DP)    , intent(in ) :: gp    (0:imax-1, 1:jmax, 1:kmax)
    real(DP)    , intent(in ) :: gt    (0:imax-1, 1:jmax, 1:kmax)
    real(DP)    , intent(in ) :: gts   (0:imax-1, 1:jmax)
    real(DP)    , intent(in ) :: dod067(0:imax-1, 1:jmax, 0:kmax)
    real(DP)    , intent(in ) :: qerat
    real(DP)    , intent(in ) :: ssa
    real(DP)    , intent(in ) :: emis  (0:imax-1, 1:jmax)

    real(DP)    , intent(out) :: gr15mnetflh    (0:imax-1, 1:jmax, 0:kmax)
    real(DP)    , intent(out) :: gdr15mnetfldtsh(0:imax-1, 1:jmax, 0:kmax, 0:1)

    real(DP)    , intent(out) :: gor (0:imax-1, 1:jmax), goru (0:imax-1, 1:jmax), gord (0:imax-1, 1:jmax), gsr (0:imax-1, 1:jmax), gsru (0:imax-1, 1:jmax), gsrd (0:imax-1, 1:jmax)


    !
    ! local variables
    !
    real(DP) :: gth(0:imax-1, 1:jmax, 0:kmax)
    real(DP) :: mmmass (0:imax-1, 1:jmax, 1:kmax)
    real(DP) :: gvmr   (0:imax-1, 1:jmax, 1:kmax, 1:nras+nrps )
    real(DP) :: ac     (0:imax-1, 1:jmax, 1:kmax, 1:nras      )
    real(DP) :: pfh    (0:imax-1, 1:jmax, 0:kmax)
    real(DP) :: pfs    (0:imax-1, 1:jmax)

    real(DP) :: pfs_for_gradcalc(0:imax-1, 1:jmax)

    real(DP) :: weight_integral
    integer  :: ig, iband

    integer  :: i, j, k, l, m, n
    integer  :: k2

    !
    ! dod      : dust optical depth
    !
    real(DP) :: gdod(0:imax-1, 1:jmax, 0:kmax)

    !
    ! local variables for pfint
    !
    integer     , parameter :: divnum = 3
    real(DP)    , parameter :: wn1 = 500.0d2, wn2 = 850.0d2

    real(DP) :: minp, maxp

    integer :: iband_reserve
    real(DP)     :: glnps(0:imax-1, 1:jmax)
    real(DP)     :: glnp (0:imax-1, 1:jmax, 1:kmax )
    integer :: jj    (0:imax-1, 1:jmax, 1:kmax), kk    (0:imax-1, 1:jmax, 1:kmax ), jjs   (0:imax-1, 1:jmax)        , kks   (0:imax-1, 1:jmax)


    ! Surface temperature for calculation of gradient of radiative flux
    real(DP) :: gts_for_gradcalc(0:imax-1, 1:jmax)
    ! Indices for calculation of gradient of radiative flux
    integer  :: jjs_for_gradcalc(0:imax-1, 1:jmax), kks_for_gradcalc(0:imax-1, 1:jmax)

    real(DP)     :: pfrh(0:imax-1, 1:jmax, 0:kmax)
    real(DP)     :: pfr (0:imax-1, 1:jmax, 1:kmax)
    real(DP)     :: pfrs(0:imax-1, 1:jmax)

    logical, save :: FlagCalcTrans

    data FlagCalcTrans / .false. /


    k = 0
    do j = 1, jmax
      do i = 0, imax-1
!!$        gth(i,j,k) = gt(i,j,k+1)
        gth(i,j,k) = ( gt(i,j,2) - gt(i,j,1) ) / log( gp (i,j,2) / gp(i,j,1) ) * log( gph(i,j,k) / gp(i,j,1) ) + gt(i,j,1)
      end do
    end do
    do k = 1, kmax-1
      do j = 1, jmax
        do i = 0, imax-1
          gth(i,j,k) = ( gt(i,j,k+1) - gt(i,j,k) ) / log( gp (i,j,k+1) / gp(i,j,k) ) * log( gph(i,j,k  ) / gp(i,j,k) ) + gt(i,j,k)
        end do
      end do
    end do
    k = kmax
    do j = 1, jmax
      do i = 0, imax-1
        gth(i,j,k) = gt(i,j,k)
      end do
    end do


!!$    do k = 1, km*nvr+1
!!$      do ij = ijs, ije
!!$        gph_f( ij, 1, k ) = sgmh_f( k ) * gph( ij, 1, km+1 )
!!$      end do
!!$    end do
!!$    call calc_lnp( im, jm, km*nvr+1, gph_f , glnph_f , ijs, ije )
!!$
!!$    call increase_vreso_boundary_arr3d( im, jm, km, nvr, gth, gth_f, &
!!$      & "linear", ijs, ije )



    !
    ! Calculation of transmission
    !
    if( ( .not. FlagCalcTrans ) .or. ( time - dble( int( time / rad15mint )  ) * rad15mint ) .lt. dt ) then

      FlagCalcTrans = .true.

!!$      call MessageNotify( 'M', module_name, 'Transmission is calculated.' )

      !
      ! Calculation of "absorption" dust optical depth
      ! This formulation is obtained from Forget et al. [1999].
      !
      do k = 0, kmax
        do j = 1, jmax
          do i = 0, imax-1
            gdod(i,j,k) = ( 1.0d0 - ssa ) * dod067(i,j,k) * qerat
          end do
        end do
      end do

!!$      call increase_vreso_boundary_arr3d( im, jm, km, nvr, gdod, gdod_f, &
!!$        & "log", ijs, ije )


      !
      ! check pressure
      !
      minp = 1.0d100
      maxp = 0.0d0
      do j = 1, jmax
        do i = 0, imax-1
          minp = min( minp, gp(i,j,kmax) )
          maxp = max( maxp, gp(i,j,1   ) )
        end do
      end do
      if( ckdp(1)%lnp(1) .gt. log(minp) ) then
        write( 6, * ) 'MARS: pressure is too small.'
        write( 6, * ) minp, exp(ckdp(1)%lnp(1))
        stop
      end if
      if( ckdp(1)%lnp(ckdp(1)%nlnp) .lt. log(maxp) ) then
        write( 6, * ) 'MARS: pressure is too large.'
        write( 6, * ) maxp, exp(ckdp(1)%lnp(ckdp(1)%nlnp))
        stop
      end if


      do k = 1, kmax
        do j = 1, jmax
          do i = 0, imax-1
            mmmass(i,j,k) = 43.5d0 * amu
          end do
        end do
      end do
      do n = 1, nras + nrps
        do k = 1, kmax
          do j = 1, jmax
            do i = 0, imax-1
              gvmr(i,j,k,n) = vmr_co2
            end do
          end do
        end do
      end do


!!$      do n = 1, nras + nrps
!!$        call increase_vreso_b2m_arr3d( im, jm, km, nvr, &
!!$          & gvmrh(:,:,:,n), gvmr_f(:,:,:,n), "log", &
!!$          & ijs, ije )
!!$      end do
!!$      call increase_vreso_b2m_arr3d( im, jm, km, nvr, &
!!$        & mmmassh, mmmass_f, "linear", &
!!$        & ijs, ije )

!!$      call calc_lnp( im, jm, km+1    , gph   , glnph   , ijs, ije )
      call calc_lnp( gp, glnp )
      glnps(:,:) = log( gph(:,:,0) )


      !
      ! initialization
      !
      do k = 0, kmax
        do j = 1, jmax
          do i = 0, imax-1
            trans_i2i_toa(i,j,k) = 0.0d0         ! f_{1/2}    T_{k+1/2,1/2}
            trans_i2i_boa(i,j,k) = 0.0d0         ! f_{km+1/2} T_{k+1/2,km+1/2}
            trans_i2i_s  (i,j,k) = 0.0d0         ! f_{s}      T_{k+1/2,km+1/2}
          end do
        end do
      end do
      do k2 = 1, kmax
        do k = 0, kmax
          do j = 1, jmax
            do i = 0, imax-1
              trans_i2m_uli(i,j,k,k2) = 0.0d0
              trans_i2m_lli(i,j,k,k2) = 0.0d0
            end do
          end do
        end do
      end do


      !
      ! loop for wavenumber
      !

      iband_reserve = 0

      do m = 1, nwnl

        call m2ckdpindices( m, ig, iband )


        if( iband .ne. iband_reserve ) then
          call findindices3D( gt, glnp, iband, jj, kk )
          call findindices2D( gts, glnps, iband, jjs, kks )

          iband_reserve = iband
        end if


        ! IMPORTANT!
        ! This loop for n is confusing.
        ! We have to reconsider about it.
        ! Maybe, the component of ckdp structure has to be reconsidered. 
        ! Now, it cannot include multiple radiatively active species. 
        ! (yot, 2010/09/12)
        !
        do n = 1, nras
          call getlnac_givenindices( gt, glnp, jj, kk, ig, iband, ac(:,:,:,n) )
        end do
        do n = 1, nras
          do k = 1, kmax
            do j = 1, jmax
              do i = 0, imax-1
                ac(i,j,k,n) = exp( ac(i,j,k,n) )
              end do
            end do
          end do
        end do

!!$        do n = 1, nras
!!$          call increase_vreso_b2m_arr3d( im, jm, km, nvr, &
!!$            & ach(:,:,:,n), ac_f(:,:,:,n), "log", &
!!$            & ijs, ije )
!!$        end do


        call calc_trans_mp_arr3d( nras, nrps, gph, gvmr, mmmass, ac, gdod, trans )


        call getpfr_givenindices3D( gt, glnp, jj, kk, ig, iband, pfr )
        pfrh(:,:,0) = pfr(:,:,1)
        do k = 1, kmax-1
          pfrh(:,:,k) = ( pfr(:,:,k) + pfr(:,:,k+1) ) * 0.5_DP
        end do
        pfrh(:,:,kmax) = pfr(:,:,kmax)
        call getpfr_givenindices2D( gts, glnps, jjs, kks, ig, iband, pfrs )

        do k = 0, kmax
          do j = 1, jmax
            do i = 0, imax-1
              trans_i2i_toa(i,j,k) = trans_i2i_toa(i,j,k) + trans(i,j,k,kmax) * pfrh(i,j,kmax) * ckdp(iband)%weight(ig)
              trans_i2i_boa(i,j,k) = trans_i2i_boa(i,j,k) + trans(i,j,k,0) * pfrh(i,j,0) * ckdp(iband)%weight(ig)
              trans_i2i_s  (i,j,k) = trans_i2i_s  (i,j,k) + trans(i,j,k,0) * pfrs(i,j) * ckdp(iband)%weight(ig)
            end do
          end do
        end do

        do k2 = 1, kmax
          do k = 0, kmax
            do j = 1, jmax
              do i = 0, imax-1
                trans_i2m_uli(i,j,k,k2) = trans_i2m_uli(i,j,k,k2) + ( trans(i,j,k,k2-1) + trans(i,j,k,k2) ) * 0.5d0 * pfrh(i,j,k2  ) * ckdp(iband)%weight(ig)
                trans_i2m_lli(i,j,k,k2) = trans_i2m_lli(i,j,k,k2) + ( trans(i,j,k,k2-1) + trans(i,j,k,k2) ) * 0.5d0 * pfrh(i,j,k2-1) * ckdp(iband)%weight(ig)
              end do
            end do
          end do
        end do


      end do


!!$      call rad15m_rv_put_newscheme2006( time, gt, gph, gp, gts, gdod, ijs, ije )

    else

      if ( trans_i2i_toa(1,1,1) > 1.0d99 ) then
        write( 6, * ) 'transmission function would not be calculated.'
        stop
      end if

    end if


    ! Is this OK?
    iband = 1

    call getpf_arr3d_norat( gth, gts, iband, pfh, pfs )


    call calc_rteq_use_meantrans_arr3d( ( ckdp(iband)%wnbnds(2) - ckdp(iband)%wnbnds(1) ), emis, trans_i2i_toa, trans_i2i_boa, trans_i2i_s, trans_i2m_lli, trans_i2m_uli, pfh, pfs, gr15mnetflh )


    do l = 0, 1
      do k = 0, kmax
        do j = 1, jmax
          do i = 0, imax-1
            gdr15mnetfldtsh(i,j,k,l) = 0.0_DP
          end do
        end do
      end do
    end do


!!$    do k = kmax, 0, -1
!!$      write( 6, * ) gph(0,1,k), gr15mnetflh(0,1,k)
!!$    end do
!!$    stop


!!$      ij = ( ije - ijs + 1 ) / 2
!!$      k  = km + 1
!!$!      write( 6, * ) 'RAD15M : ', gr15mnetflh(ij,1,k), rad_gdr15mnetfldtsh0(ij,1,k-1) * ( gts(ij,1) - rad_gtsbase(ij,1,1) ), emis(ij,1) * 5.67d-8 * gts(ij,1)**4, gts(ij,1)
!!$!      write( 61, * ) gr15mnetflh(ij,1,k-1), rad_gdr15mnetfldtsh0(ij,1,k-1) * ( gts(ij,1) - rad_gtsbase(ij,1,1) ), emis(ij,1) * 5.67d-8 * gts(ij,1)**4, gts(ij,1), gt(ij,1,km), gt(ij,1,km-1), gt(ij,1,km-2), gt(ij,1,km-3), gt(ij,1,km-4), gt(ij,1,km-5)
!!$      write( 6, * ) 'RAD15M : ', gr15mnetflh(ij,1,k), &
!!$        & emis(ij,1) * 5.67d-8 * gts(ij,1)**4, gts(ij,1)
!!$      write( 61, * ) gr15mnetflh(ij,1,k-1), &
!!$        & emis(ij,1) * 5.67d-8 * gts(ij,1)**4, gts(ij,1), &
!!$        & gt(ij,1,km), gt(ij,1,km-1), gt(ij,1,km-2), gt(ij,1,km-3), &
!!$        & gt(ij,1,km-4), gt(ij,1,km-5)
!!$      call flush( 61 )

    !
    ! output variables
    !
!!$    do j = 1, jmax
!!$      do i = 0, imax-1
!!$        goru(i,j) = uwflh_sum(i,j,kmax)
!!$        gord(i,j) = 0.0d0
!!$        gsru(i,j) = uwflh_sum(i,j,0)
!!$        gsrd(i,j) = dwflh_sum(i,j,0)
!!$        gor (i,j) = goru(i,j) - gord(i,j)
!!$        gsr (i,j) = gsru(i,j) - gsrd(i,j)
!!$      end do
!!$    end do



  end subroutine rad15m_lowatm_newscheme2006
rad15mint
Variable :
rad15mint :real(DP) , save
rad_gdod
Variable :
rad_gdod( :, :, : ) :real(DP), allocatable, save
rad_gp
Variable :
rad_gp( :, :, : ) :real(DP), allocatable, save
rad_gph
Variable :
rad_gph( :, :, : ) :real(DP), allocatable, save
rad_gt
Variable :
rad_gt( :, :, : ) :real(DP), allocatable, save
rad_gts
Variable :
rad_gts( :, :, : ) :real(DP), allocatable, save
rad_time
Variable :
rad_time :real(DP), save
sw_prep_rv
Variable :
sw_prep_rv :logical , save
trans
Variable :
trans(:,:,:,:) :real(DP) , allocatable, save
:
!$ real(DP) , allocatable, save :gdod_f ( :, :, : )

!$

!$ real(DP) , allocatable, save :uwflh_f( :, :, : ), dwflh_f( :, :, : )
trans_i2i_boa
Variable :
trans_i2i_boa(:,:,:) :real(DP) , allocatable, save
: lower layer interface
trans_i2i_s
Variable :
trans_i2i_s(:,:,:) :real(DP) , allocatable, save
: lower layer interface
trans_i2i_toa
Variable :
trans_i2i_toa(:,:,:) :real(DP) , allocatable, save
: lower layer interface
trans_i2m_lli
Variable :
trans_i2m_lli(:,:,:,:) :real(DP) , allocatable, save
: lower layer interface
trans_i2m_uli
Variable :
trans_i2m_uli(:,:,:,:) :real(DP) , allocatable, save
: lower layer interface
trans_res
Variable :
trans_res(:,:,:,:) :real(DP), allocatable, save
version
Constant :
version = ’$Name: dcpam5-20110225-4 $’ // ’$Id: radiation_dcpam_M_15m.f90,v 1.1 2010-09-18 01:45:09 yot Exp $’ :character(*), parameter
: モジュールのバージョン Module version
vmr_co2
Constant :
vmr_co2 = 0.95d0 :real(DP), parameter