SCALE-RM
scale_atmos_phy_rd_common.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
14 !-------------------------------------------------------------------------------
16  !-----------------------------------------------------------------------------
17  !
18  !++ used modules
19  !
20  use scale_precision
21  use scale_stdio
22  use scale_prof
24  use scale_tracer
25  !-----------------------------------------------------------------------------
26  implicit none
27  private
28  !-----------------------------------------------------------------------------
29  !
30  !++ Public procedure
31  !
32  public :: atmos_phy_rd_heating
33 
34  !-----------------------------------------------------------------------------
35  !
36  !++ Public parameters & variables
37  !
38  ! for direction
39  integer, public, parameter :: i_up = 1
40  integer, public, parameter :: i_dn = 2
41  ! for band region
42  integer, public, parameter :: i_lw = 1
43  integer, public, parameter :: i_sw = 2
44  ! for band region
45  integer, public, parameter :: i_direct = 1
46  integer, public, parameter :: i_diffuse = 2
47 
48  !-----------------------------------------------------------------------------
49  !
50  !++ Private procedure
51  !
52  !-----------------------------------------------------------------------------
53  !
54  !++ Private parameters & variables
55  !
56  !-----------------------------------------------------------------------------
57 contains
58  !-----------------------------------------------------------------------------
60  subroutine atmos_phy_rd_heating( &
61  flux_rad, &
62  DENS, &
63  RHOT, &
64  QTRC, &
65  FZ, &
66  dt, &
67  TEMP_t, &
68  RHOT_t )
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
148  end subroutine atmos_phy_rd_heating
149 
150 end module scale_atmos_phy_rd_common
integer, public is
start point of inner domain: x, local
subroutine, public atmos_phy_rd_heating(flux_rad, DENS, RHOT, QTRC, FZ, dt, TEMP_t, RHOT_t)
Calc heating rate.
integer, public je
end point of inner domain: y, local
integer, parameter, public i_direct
module STDIO
Definition: scale_stdio.F90:12
integer, public ke
end point of inner domain: z, local
integer, public qa
integer, parameter, public i_diffuse
integer, parameter, public i_lw
module grid index
integer, parameter, public i_sw
module TRACER
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
module profiler
Definition: scale_prof.F90:10
integer, public ie
end point of inner domain: x, local
module ATMOSPHERE / Thermodynamics
module PRECISION
module ATMOSPHERE / Physics Radiation
integer, parameter, public i_up
integer, public ja
of y whole cells (local, with HALO)