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 70 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, scale_grid_index::ks, scale_tracer::tracer_cv, scale_tracer::tracer_mass, and scale_tracer::tracer_r.

Referenced by mod_atmos_phy_rd_driver::atmos_phy_rd_driver().

70  use scale_tracer, only: &
71  tracer_r, &
72  tracer_cv, &
74  use scale_atmos_thermodyn, only: &
75  thermodyn_qd => atmos_thermodyn_qd, &
76  thermodyn_cv => atmos_thermodyn_cv, &
77  thermodyn_rhoe => atmos_thermodyn_rhoe, &
78  thermodyn_rhot => atmos_thermodyn_rhot
79  implicit none
80 
81  real(RP), intent(in) :: flux_rad(KA,IA,JA,2,2)
82  real(RP), intent(in) :: DENS (KA,IA,JA)
83  real(RP), intent(in) :: RHOT (KA,IA,JA)
84  real(RP), intent(in) :: QTRC (KA,IA,JA,QA)
85  real(RP), intent(in) :: FZ (0:KA,IA,JA)
86  real(DP), intent(in) :: dt
87  real(RP), intent(out) :: TEMP_t (KA,IA,JA,3)
88  real(RP), intent(out) :: RHOT_t (KA,IA,JA)
89 
90  real(RP) :: RHOE (KA,IA,JA)
91  real(RP) :: RHOE_t(KA,IA,JA,2)
92  real(RP) :: RHOT1 (KA,IA,JA)
93  real(RP) :: QDRY (KA,IA,JA)
94  real(RP) :: CVtot (KA,IA,JA)
95 
96  integer :: k, i, j
97  !---------------------------------------------------------------------------
98 
99  call thermodyn_rhoe( rhoe(:,:,:), & ! [OUT]
100  rhot(:,:,:), & ! [IN]
101  qtrc(:,:,:,:), & ! [IN]
102  tracer_cv(:), & ! [IN]
103  tracer_r(:), & ! [IN]
104  tracer_mass(:) ) ! [IN]
105 
106  !$omp parallel do default(none) &
107  !$omp shared(JS,JE,IS,IE,KS,KE,RHOE_t,flux_rad,FZ,dt,RHOE) &
108  !$omp private(i,j,k) OMP_SCHEDULE_ collapse(2)
109  do j = js, je
110  do i = is, ie
111  do k = ks, ke
112 
113  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) ) &
114  - ( flux_rad(k,i,j,i_lw,i_up) - flux_rad(k-1,i,j,i_lw,i_up) ) &
115  ) / ( fz(k,i,j) - fz(k-1,i,j) )
116 
117  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) ) &
118  - ( flux_rad(k,i,j,i_sw,i_up) - flux_rad(k-1,i,j,i_sw,i_up) ) &
119  ) / ( fz(k,i,j) - fz(k-1,i,j) )
120 
121  rhoe(k,i,j) = rhoe(k,i,j) + dt * ( rhoe_t(k,i,j,i_lw) + rhoe_t(k,i,j,i_sw) )
122 
123  enddo
124  enddo
125  enddo
126 
127  call thermodyn_rhot( rhot1(:,:,:), & ! [OUT]
128  rhoe(:,:,:), & ! [IN]
129  qtrc(:,:,:,:), & ! [IN]
130  tracer_cv(:), & ! [IN]
131  tracer_r(:), & ! [IN]
132  tracer_mass(:) ) ! [IN]
133 
134  ! update rhot
135  do j = js, je
136  do i = is, ie
137  do k = ks, ke
138  rhot_t(k,i,j) = ( rhot1(k,i,j) - rhot(k,i,j) ) / dt
139  enddo
140  enddo
141  enddo
142 
143  call thermodyn_qd( qdry(:,:,:), & ! [OUT]
144  qtrc(:,:,:,:), & ! [IN]
145  tracer_mass(:) ) ! [IN]
146 
147  call thermodyn_cv( cvtot(:,:,:), & ! [OUT]
148  qtrc(:,:,:,:), & ! [IN]
149  tracer_cv(:), & ! [IN]
150  qdry(:,:,:) ) ! [IN]
151 
152  do j = js, je
153  do i = is, ie
154  do k = ks, ke
155  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]
156  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]
157 
158  temp_t(k,i,j,3) = temp_t(k,i,j,i_lw) + temp_t(k,i,j,i_sw)
159  enddo
160  enddo
161  enddo
162 
163  return
integer, public is
start point of inner domain: x, local
integer, public je
end point of inner domain: y, local
real(rp), dimension(qa_max), public tracer_r
integer, public ke
end point of inner domain: z, local
real(rp), dimension(qa_max), public tracer_cv
integer, parameter, public i_lw
integer, parameter, public i_sw
module TRACER
integer, parameter, public i_dn
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
real(rp), dimension(qa_max), public tracer_mass
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 46 of file scale_atmos_phy_rd_common.F90.

Referenced by mod_atmos_phy_rd_driver::atmos_phy_rd_driver(), and scale_atmos_phy_rd_offline::atmos_phy_rd_offline().

46  integer, public, parameter :: I_direct = 1

◆ i_diffuse

integer, parameter, public scale_atmos_phy_rd_common::i_diffuse = 2

Definition at line 47 of file scale_atmos_phy_rd_common.F90.

Referenced by mod_atmos_phy_rd_driver::atmos_phy_rd_driver(), and scale_atmos_phy_rd_offline::atmos_phy_rd_offline().

47  integer, public, parameter :: I_diffuse = 2