SCALE-RM
Functions/Subroutines | Variables
scale_atmos_phy_rd_common Module Reference

module ATMOSPHERE / Physics Radiation More...

Functions/Subroutines

subroutine, public atmos_phy_rd_heating (flux_rad, DENS, RHOT, QTRC, FZ, dt, TEMP_t, RHOT_t)
 Calc heating rate. More...
 

Variables

integer, parameter, public i_up = 1
 
integer, parameter, public i_dn = 2
 
integer, parameter, public i_lw = 1
 
integer, parameter, public i_sw = 2
 
integer, parameter, public i_direct = 1
 
integer, parameter, public i_diffuse = 2
 

Detailed Description

module ATMOSPHERE / Physics Radiation

Description
Common module for Radiation Heating Rate
Author
Team SCALE
History
  • 2013-02-06 (H.Yashiro) [new]

Function/Subroutine Documentation

◆ atmos_phy_rd_heating()

subroutine, public scale_atmos_phy_rd_common::atmos_phy_rd_heating ( real(rp), dimension(ka,ia,ja,2,2), intent(in)  flux_rad,
real(rp), dimension (ka,ia,ja), intent(in)  DENS,
real(rp), dimension (ka,ia,ja), intent(in)  RHOT,
real(rp), dimension (ka,ia,ja,qa), intent(in)  QTRC,
real(rp), dimension (0:ka,ia,ja), intent(in)  FZ,
real(dp), intent(in)  dt,
real(rp), dimension (ka,ia,ja,3), intent(out)  TEMP_t,
real(rp), dimension (ka,ia,ja), intent(out)  RHOT_t 
)

Calc heating rate.

Definition at line 69 of file scale_atmos_phy_rd_common.F90.

References i_dn, i_lw, i_sw, i_up, scale_grid_index::ie, scale_grid_index::is, scale_grid_index::je, scale_grid_index::js, scale_grid_index::ke, and scale_grid_index::ks.

Referenced by mod_atmos_phy_rd_driver::atmos_phy_rd_driver().

69  use scale_atmos_thermodyn, only: &
70  thermodyn_qd => atmos_thermodyn_qd, &
71  thermodyn_cv => atmos_thermodyn_cv, &
72  thermodyn_rhoe => atmos_thermodyn_rhoe, &
73  thermodyn_rhot => atmos_thermodyn_rhot
74  implicit none
75 
76  real(RP), intent(in) :: flux_rad(ka,ia,ja,2,2)
77  real(RP), intent(in) :: dens (ka,ia,ja)
78  real(RP), intent(in) :: rhot (ka,ia,ja)
79  real(RP), intent(in) :: qtrc (ka,ia,ja,qa)
80  real(RP), intent(in) :: fz (0:ka,ia,ja)
81  real(DP), intent(in) :: dt
82  real(RP), intent(out) :: temp_t (ka,ia,ja,3)
83  real(RP), intent(out) :: rhot_t (ka,ia,ja)
84 
85  real(RP) :: rhoe (ka,ia,ja)
86  real(RP) :: rhoe_t(ka,ia,ja,2)
87  real(RP) :: rhot1 (ka,ia,ja)
88  real(RP) :: qdry (ka,ia,ja)
89  real(RP) :: cvtot (ka,ia,ja)
90 
91  integer :: k, i, j
92  !---------------------------------------------------------------------------
93 
94  call thermodyn_rhoe( rhoe(:,:,:), & ! [OUT]
95  rhot(:,:,:), & ! [IN]
96  qtrc(:,:,:,:) ) ! [IN]
97 
98  do j = js, je
99  do i = is, ie
100  do k = ks, ke
101 
102  rhoe_t(k,i,j,i_lw) = ( ( flux_rad(k,i,j,i_lw,i_dn) - flux_rad(k-1,i,j,i_lw,i_dn) ) &
103  - ( flux_rad(k,i,j,i_lw,i_up) - flux_rad(k-1,i,j,i_lw,i_up) ) &
104  ) / ( fz(k,i,j) - fz(k-1,i,j) )
105 
106  rhoe_t(k,i,j,i_sw) = ( ( flux_rad(k,i,j,i_sw,i_dn) - flux_rad(k-1,i,j,i_sw,i_dn) ) &
107  - ( flux_rad(k,i,j,i_sw,i_up) - flux_rad(k-1,i,j,i_sw,i_up) ) &
108  ) / ( fz(k,i,j) - fz(k-1,i,j) )
109 
110  rhoe(k,i,j) = rhoe(k,i,j) + dt * ( rhoe_t(k,i,j,i_lw) + rhoe_t(k,i,j,i_sw) )
111 
112  enddo
113  enddo
114  enddo
115 
116  call thermodyn_rhot( rhot1(:,:,:), & ! [OUT]
117  rhoe(:,:,:), & ! [IN]
118  qtrc(:,:,:,:) ) ! [IN]
119 
120  ! update rhot
121  do j = js, je
122  do i = is, ie
123  do k = ks, ke
124  rhot_t(k,i,j) = ( rhot1(k,i,j) - rhot(k,i,j) ) / dt
125  enddo
126  enddo
127  enddo
128 
129  call thermodyn_qd( qdry(:,:,:), & ! [OUT]
130  qtrc(:,:,:,:) ) ! [IN]
131 
132  call thermodyn_cv( cvtot(:,:,:), & ! [OUT]
133  qtrc(:,:,:,:), & ! [IN]
134  qdry(:,:,:) ) ! [IN]
135 
136  do j = js, je
137  do i = is, ie
138  do k = ks, ke
139  temp_t(k,i,j,i_lw) = rhoe_t(k,i,j,i_lw) / dens(k,i,j) / cvtot(k,i,j) * 86400.0_rp ! [K/day]
140  temp_t(k,i,j,i_sw) = rhoe_t(k,i,j,i_sw) / dens(k,i,j) / cvtot(k,i,j) * 86400.0_rp ! [K/day]
141 
142  temp_t(k,i,j,3) = temp_t(k,i,j,i_lw) + temp_t(k,i,j,i_sw)
143  enddo
144  enddo
145  enddo
146 
147  return
integer, public is
start point of inner domain: x, local
integer, public je
end point of inner domain: y, local
integer, public ke
end point of inner domain: z, local
integer, public qa
integer, parameter, public i_lw
integer, parameter, public i_sw
integer, public ia
of x whole cells (local, with HALO)
integer, parameter, public i_dn
integer, public ka
of z whole cells (local, with HALO)
integer, public js
start point of inner domain: y, local
integer, public ks
start point of inner domain: z, local
integer, public ie
end point of inner domain: x, local
module ATMOSPHERE / Thermodynamics
integer, parameter, public i_up
integer, public ja
of y whole cells (local, with HALO)
Here is the caller graph for this function:

Variable Documentation

◆ i_up

integer, parameter, public scale_atmos_phy_rd_common::i_up = 1

◆ i_dn

integer, parameter, public scale_atmos_phy_rd_common::i_dn = 2

◆ i_lw

integer, parameter, public scale_atmos_phy_rd_common::i_lw = 1

◆ i_sw

integer, parameter, public scale_atmos_phy_rd_common::i_sw = 2

◆ i_direct

integer, parameter, public scale_atmos_phy_rd_common::i_direct = 1

Definition at line 45 of file scale_atmos_phy_rd_common.F90.

Referenced by mod_atmos_phy_rd_driver::atmos_phy_rd_driver().

45  integer, public, parameter :: i_direct = 1
integer, parameter, public i_direct

◆ i_diffuse

integer, parameter, public scale_atmos_phy_rd_common::i_diffuse = 2

Definition at line 46 of file scale_atmos_phy_rd_common.F90.

Referenced by mod_atmos_phy_rd_driver::atmos_phy_rd_driver().

46  integer, public, parameter :: i_diffuse = 2
integer, parameter, public i_diffuse