SCALE-RM
Functions/Subroutines
scale_atmos_phy_bl_common Module Reference

module atmosphere / physics / pbl / common More...

Functions/Subroutines

subroutine, public atmos_phy_bl_tendency_tracer (KA, KS, KE, IA, IS, IE, JA, JS, JE, DENS, QTRC, SFLX_Q, Kh, MASS, CZ, FZ, F2H, DDT, TRACER_NAME, RHOQ_t)
 ATMOS_PHY_BL_tendency_tracer calculate tendency of tracers by the eddy viscosity. More...
 

Detailed Description

module atmosphere / physics / pbl / common

Description
Common routines for boundary layer turbulence
Author
Team SCALE
NAMELIST
  • No namelist group
History Output
namedescriptionunitvariable
ZFLX_{TRACER_NAME}_BL Z FLUX of DENS * {TRACER_NAME} (PBL);
{TRACER_NAME} depends on the physics schemes, e.g., QV, QC, QR.
kg/m2/s flx

Function/Subroutine Documentation

◆ atmos_phy_bl_tendency_tracer()

subroutine, public scale_atmos_phy_bl_common::atmos_phy_bl_tendency_tracer ( integer, intent(in)  KA,
integer, intent(in)  KS,
integer, intent(in)  KE,
integer, intent(in)  IA,
integer, intent(in)  IS,
integer, intent(in)  IE,
integer, intent(in)  JA,
integer, intent(in)  JS,
integer, intent(in)  JE,
real(rp), dimension (ka,ia,ja), intent(in)  DENS,
real(rp), dimension (ka,ia,ja), intent(in)  QTRC,
real(rp), dimension( ia,ja), intent(in)  SFLX_Q,
real(rp), dimension (ka,ia,ja), intent(in)  Kh,
real(rp), intent(in)  MASS,
real(rp), dimension ( ka,ia,ja), intent(in)  CZ,
real(rp), dimension (0:ka,ia,ja), intent(in)  FZ,
real(rp), dimension(ka,2,ia,ja), intent(in)  F2H,
real(dp), intent(in)  DDT,
character(len=*), intent(in)  TRACER_NAME,
real(rp), dimension(ka,ia,ja), intent(out)  RHOQ_t 
)

ATMOS_PHY_BL_tendency_tracer calculate tendency of tracers by the eddy viscosity.

Definition at line 62 of file scale_atmos_phy_bl_common.F90.

62  use scale_matrix, only: &
63  matrix_solver_tridiagonal
64  use scale_file_history, only: &
65  file_history_in
66  integer, intent(in) :: KA, KS, KE
67  integer, intent(in) :: IA, IS, IE
68  integer, intent(in) :: JA, JS, JE
69 
70  real(RP), intent(in) :: DENS (KA,IA,JA)
71  real(RP), intent(in) :: QTRC (KA,IA,JA)
72  real(RP), intent(in) :: SFLX_Q( IA,JA)
73  real(RP), intent(in) :: Kh (KA,IA,JA)
74  real(RP), intent(in) :: MASS
75  real(RP), intent(in) :: CZ ( KA,IA,JA)
76  real(RP), intent(in) :: FZ (0:KA,IA,JA)
77  real(RP), intent(in) :: F2H(KA,2,IA,JA)
78  real(DP), intent(in) :: DDT
79  character(len=*), intent(in) :: TRACER_NAME
80 
81  real(RP), intent(out) :: RHOQ_t(KA,IA,JA)
82 
83  real(RP) :: QTRC_n(KA)
84  real(RP) :: RHO (KA)
85  real(RP) :: RHOKh (KA)
86  real(RP) :: a(KA)
87  real(RP) :: b(KA)
88  real(RP) :: c(KA)
89  real(RP) :: d(KA)
90  real(RP) :: rho_h
91  real(RP) :: ap
92  real(RP) :: sf_t
93 
94  real(RP) :: flx(KA,IA,JA)
95 
96  real(RP) :: CDZ(KA)
97  real(RP) :: FDZ(KA)
98 
99  real(RP) :: dt
100 
101  integer :: KE_PBL
102  integer :: k, i, j
103 
104 #ifdef _OPENACC
105  real(RP) :: work(KS:KE,4)
106 #endif
107  dt = real( ddt, kind=rp )
108 
109  !$acc data create(flx)
110 
111 !OCL INDEPENDENT
112  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
113  !$omp shared(KA,KS,KE,IS,IE,JS,JE) &
114  !$omp shared(RHOQ_t,DENS,QTRC,SFLX_Q,Kh,MASS,CZ,FZ,F2H,DT,flx) &
115  !$omp private(QTRC_n,RHO,RHOKh,rho_h,a,b,c,d,ap,sf_t,CDZ,FDZ) &
116  !$omp private(KE_PBL,k,i,j)
117  !$acc kernels
118  !$acc loop independent collapse(2) &
119  !$acc private(QTRC_n,RHO,RHOKh,rho_h,a,b,c,d,ap,sf_t,CDZ,FDZ, &
120  !$acc KE_PBL, &
121  !$acc work)
122  do j = js, je
123  do i = is, ie
124 
125  ke_pbl = ke-1
126  !$acc loop seq
127  do k = ke-2, ks+1, -1
128  if ( kh(k,i,j) > 0.0_rp ) then
129  ke_pbl = k + 1
130  else
131  exit
132  end if
133  end do
134 
135  do k = ks, ke_pbl
136  cdz(k) = fz(k ,i,j) - fz(k-1,i,j)
137  fdz(k) = cz(k+1,i,j) - cz(k ,i,j)
138  end do
139 
140  sf_t = sflx_q(i,j) / cdz(ks)
141 
142  rho(ks) = dens(ks,i,j) + dt * sf_t * mass
143  do k = ks+1, ke_pbl
144  rho(k) = dens(k,i,j)
145  end do
146 
147  ! dens * coefficient at the half level
148  do k = ks, ke_pbl-1
149  rho_h = f2h(k,1,i,j) * dens(k+1,i,j) + f2h(k,2,i,j) * dens(k,i,j)
150  rhokh(k) = rho_h * kh(k,i,j)
151  end do
152 
153  d(ks) = ( qtrc(ks,i,j) * dens(ks,i,j) + dt * sf_t ) / rho(ks)
154  do k = ks+1, ke_pbl
155  d(k) = qtrc(k,i,j)
156  end do
157 
158  c(ks) = 0.0_rp
159  do k = ks, ke_pbl-1
160  ap = - dt * rhokh(k) / fdz(k)
161  a(k) = ap / ( rho(k) * cdz(k) )
162 #ifdef _OPENACC
163  if ( k==ks ) then
164  b(k) = - a(k) + 1.0_rp
165  else
166  b(k) = - a(k) + dt * rhokh(k-1) / ( fdz(k-1) * rho(k) * cdz(k) ) + 1.0_rp
167  end if
168 #else
169  b(k) = - a(k) - c(k) + 1.0_rp
170 #endif
171  c(k+1) = ap / ( rho(k+1) * cdz(k+1) )
172  end do
173  a(ke_pbl) = 0.0_rp
174  b(ke_pbl) = - c(ke_pbl) + 1.0_rp
175 
176  call matrix_solver_tridiagonal( &
177  ka, ks, ke_pbl, &
178 #ifdef _OPENACC
179  work(:,:), & ! (wrok)
180 #endif
181  a(:), b(:), c(:), d(:), & ! (in)
182  qtrc_n(:) ) ! (out)
183 
184  rhoq_t(ks,i,j) = ( qtrc_n(ks) * rho(ks) - qtrc(ks,i,j) * dens(ks,i,j) ) / dt - sf_t
185  do k = ks+1, ke_pbl
186  rhoq_t(k,i,j) = ( qtrc_n(k) - qtrc(k,i,j) ) * rho(k) / dt
187  end do
188  do k = ke_pbl+1, ke
189  rhoq_t(k,i,j) = 0.0_rp
190  end do
191 
192  do k = ks, ke_pbl-1
193  flx(k,i,j) = - rhokh(k) * ( qtrc_n(k+1) - qtrc_n(k) ) / fdz(k)
194  end do
195 
196  end do
197  end do
198  !$acc end kernels
199 
200  call file_history_in(flx(:,:,:), 'ZFLX_'//trim(tracer_name)//'_BL', 'Z FLUX of DENS * '//trim(tracer_name)//' (PBL)', 'kg/m2/s', fill_halo=.true.)
201 
202  !$acc end data
203 
204  return

References scale_atmos_grid_cartesc_index::ie, scale_atmos_grid_cartesc_index::is, scale_atmos_grid_cartesc_index::je, scale_atmos_grid_cartesc_index::js, scale_tracer::k, scale_atmos_grid_cartesc_index::ka, scale_atmos_grid_cartesc_index::ke, scale_atmos_grid_cartesc_index::ks, scale_tracer::mass, scale_precision::rp, and scale_tracer::tracer_name.

Referenced by mod_atmos_phy_bl_driver::atmos_phy_bl_driver_calc_tendency().

Here is the caller graph for this function:
scale_file_history
module file_history
Definition: scale_file_history.F90:15
scale_matrix
module MATRIX
Definition: scale_matrix.F90:17