SCALE-RM
scale_atmos_phy_rd_common.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
10 !-------------------------------------------------------------------------------
11 #include "scalelib.h"
13  !-----------------------------------------------------------------------------
14  !
15  !++ used modules
16  !
17  use scale_precision
18  use scale_io
19  use scale_prof
20  !-----------------------------------------------------------------------------
21  implicit none
22  private
23  !-----------------------------------------------------------------------------
24  !
25  !++ Public procedure
26  !
28 
29  !-----------------------------------------------------------------------------
30  !
31  !++ Public parameters & variables
32  !
33  ! for direction
34  integer, public, parameter :: i_up = 1
35  integer, public, parameter :: i_dn = 2
36  ! for band region
37  integer, public, parameter :: i_lw = 1
38  integer, public, parameter :: i_sw = 2
39 
40  !-----------------------------------------------------------------------------
41  !
42  !++ Private procedure
43  !
44  !-----------------------------------------------------------------------------
45  !
46  !++ Private parameters & variables
47  !
48  !-----------------------------------------------------------------------------
49 contains
50  !-----------------------------------------------------------------------------
52  subroutine atmos_phy_rd_calc_heating( &
53  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
54  flux_rad, &
55  DENS, TEMP, &
56  CVtot, &
57  FZ, &
58  RHOH, &
59  TEMP_t )
60  implicit none
61  integer, intent(in) :: ka, ks, ke
62  integer, intent(in) :: ia, is, ie
63  integer, intent(in) :: ja, js, je
64 
65  real(rp), intent(in) :: flux_rad(ka,ia,ja,2,2)
66  real(rp), intent(in) :: dens (ka,ia,ja)
67  real(rp), intent(in) :: temp (ka,ia,ja)
68  real(rp), intent(in) :: cvtot (ka,ia,ja)
69  real(rp), intent(in) :: fz (0:ka,ia,ja)
70 
71  real(rp), intent(out) :: rhoh(ka,ia,ja)
72 
73  real(rp), intent(out), optional :: temp_t(ka,ia,ja,3)
74 
75  real(rp) :: rhoh_lw, rhoh_sw
76 
77  integer :: k, i, j
78  !---------------------------------------------------------------------------
79 
80  if ( present(temp_t) ) then
81 
82  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
83  !$omp private(i,j,k, &
84  !$omp RHOH_LW,RHOH_SW) &
85  !$omp shared(RHOH,TEMP_t,flux_rad,DENS,CVtot,FZ, &
86  !$omp KS,KE,IS,IE,JS,JE)
87  !$acc kernels
88  do j = js, je
89  do i = is, ie
90  do k = ks, ke
91 
92  rhoh_lw = ( ( flux_rad(k,i,j,i_lw,i_dn) - flux_rad(k-1,i,j,i_lw,i_dn) ) &
93  - ( flux_rad(k,i,j,i_lw,i_up) - flux_rad(k-1,i,j,i_lw,i_up) ) &
94  ) / ( fz(k,i,j) - fz(k-1,i,j) )
95 
96  rhoh_sw = ( ( flux_rad(k,i,j,i_sw,i_dn) - flux_rad(k-1,i,j,i_sw,i_dn) ) &
97  - ( flux_rad(k,i,j,i_sw,i_up) - flux_rad(k-1,i,j,i_sw,i_up) ) &
98  ) / ( fz(k,i,j) - fz(k-1,i,j) )
99 
100  rhoh(k,i,j) = rhoh_lw + rhoh_sw
101 
102  temp_t(k,i,j,i_lw) = rhoh_lw / dens(k,i,j) / cvtot(k,i,j) * 86400.0_rp ! [K/day]
103  temp_t(k,i,j,i_sw) = rhoh_sw / dens(k,i,j) / cvtot(k,i,j) * 86400.0_rp ! [K/day]
104 
105  temp_t(k,i,j,3) = temp_t(k,i,j,i_lw) + temp_t(k,i,j,i_sw)
106 
107  enddo
108  enddo
109  enddo
110  !$acc end kernels
111 
112  else
113 
114  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
115  !$omp private(i,j,k, &
116  !$omp RHOH_LW,RHOH_SW) &
117  !$omp shared(RHOH,flux_rad,FZ, &
118  !$omp KS,KE,IS,IE,JS,JE)
119  !$acc kernels
120  do j = js, je
121  do i = is, ie
122  do k = ks, ke
123 
124  rhoh_lw = ( ( flux_rad(k,i,j,i_lw,i_dn) - flux_rad(k-1,i,j,i_lw,i_dn) ) &
125  - ( flux_rad(k,i,j,i_lw,i_up) - flux_rad(k-1,i,j,i_lw,i_up) ) &
126  ) / ( fz(k,i,j) - fz(k-1,i,j) )
127 
128  rhoh_sw = ( ( flux_rad(k,i,j,i_sw,i_dn) - flux_rad(k-1,i,j,i_sw,i_dn) ) &
129  - ( flux_rad(k,i,j,i_sw,i_up) - flux_rad(k-1,i,j,i_sw,i_up) ) &
130  ) / ( fz(k,i,j) - fz(k-1,i,j) )
131 
132  rhoh(k,i,j) = rhoh_lw + rhoh_sw
133 
134  enddo
135  enddo
136  enddo
137  !$acc end kernels
138 
139  end if
140 
141  return
142  end subroutine atmos_phy_rd_calc_heating
143 
144 end module scale_atmos_phy_rd_common
scale_atmos_phy_rd_common::atmos_phy_rd_calc_heating
subroutine, public atmos_phy_rd_calc_heating(KA, KS, KE, IA, IS, IE, JA, JS, JE, flux_rad, DENS, TEMP, CVtot, FZ, RHOH, TEMP_t)
Calc heating rate.
Definition: scale_atmos_phy_rd_common.F90:60
scale_precision
module PRECISION
Definition: scale_precision.F90:14
scale_atmos_phy_rd_common::i_dn
integer, parameter, public i_dn
Definition: scale_atmos_phy_rd_common.F90:35
scale_atmos_phy_rd_common::i_lw
integer, parameter, public i_lw
Definition: scale_atmos_phy_rd_common.F90:37
scale_atmos_phy_rd_common::i_up
integer, parameter, public i_up
Definition: scale_atmos_phy_rd_common.F90:34
scale_precision::rp
integer, parameter, public rp
Definition: scale_precision.F90:41
scale_io
module STDIO
Definition: scale_io.F90:10
scale_prof
module profiler
Definition: scale_prof.F90:11
scale_atmos_phy_rd_common::i_sw
integer, parameter, public i_sw
Definition: scale_atmos_phy_rd_common.F90:38
scale_atmos_phy_rd_common
module atmosphere / physics / radiation / common
Definition: scale_atmos_phy_rd_common.F90:12