SCALE-RM
scale_atmos_phy_rd_common.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
14 !-------------------------------------------------------------------------------
15 #include "inc_openmp.h"
17  !-----------------------------------------------------------------------------
18  !
19  !++ used modules
20  !
21  use scale_precision
22  use scale_stdio
23  use scale_prof
25  use scale_tracer
26  !-----------------------------------------------------------------------------
27  implicit none
28  private
29  !-----------------------------------------------------------------------------
30  !
31  !++ Public procedure
32  !
33  public :: atmos_phy_rd_heating
34 
35  !-----------------------------------------------------------------------------
36  !
37  !++ Public parameters & variables
38  !
39  ! for direction
40  integer, public, parameter :: i_up = 1
41  integer, public, parameter :: i_dn = 2
42  ! for band region
43  integer, public, parameter :: i_lw = 1
44  integer, public, parameter :: i_sw = 2
45  ! for band region
46  integer, public, parameter :: i_direct = 1
47  integer, public, parameter :: i_diffuse = 2
48 
49  !-----------------------------------------------------------------------------
50  !
51  !++ Private procedure
52  !
53  !-----------------------------------------------------------------------------
54  !
55  !++ Private parameters & variables
56  !
57  !-----------------------------------------------------------------------------
58 contains
59  !-----------------------------------------------------------------------------
61  subroutine atmos_phy_rd_heating( &
62  flux_rad, &
63  DENS, &
64  RHOT, &
65  QTRC, &
66  FZ, &
67  dt, &
68  TEMP_t, &
69  RHOT_t )
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
164  end subroutine atmos_phy_rd_heating
165 
166 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
real(rp), dimension(qa_max), public tracer_r
module STDIO
Definition: scale_stdio.F90:12
integer, public ke
end point of inner domain: z, local
integer, parameter, public i_diffuse
real(rp), dimension(qa_max), public tracer_cv
integer, parameter, public i_lw
module grid index
integer, parameter, public i_sw
module TRACER
integer, public ia
of whole cells: x, local, with HALO
integer, parameter, public i_dn
integer, public ka
of whole cells: z, 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
real(rp), dimension(qa_max), public tracer_mass
integer, public ja
of whole cells: y, local, with HALO